-
Notifications
You must be signed in to change notification settings - Fork 1
/
MATRIX_BLOCKS.f90
1148 lines (1006 loc) · 38.1 KB
/
MATRIX_BLOCKS.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
MODULE MATRIX_BLOCKS
USE SCIFOR, only: str,free_unit,zero,assert_shape,zeye,eigh,sort_quicksort
USE AUX_FUNCS
USE MATRIX_SPARSE, only:sparse_matrix,as_matrix
implicit none
private
!BLOCK MATRIX COMPONENT
type block_type
integer :: index=0
real(8),dimension(:),allocatable :: qn
real(8),dimension(:),allocatable :: E
integer,dimension(:),allocatable :: map
real(8),dimension(:,:),allocatable :: M
type(block_type),pointer :: next=>null()
end type block_type
!BLOCK MATRIX STRUCTURE
type blocks_matrix
integer :: Nblock=0
integer :: Nrow=0
integer :: Ncol=0
real(8),dimension(:),allocatable :: evalues
integer,dimension(:),allocatable :: eorder
logical :: diag=.false.
type(block_type),pointer :: root=>null()
contains
procedure,pass :: free => free_blocks_matrix
procedure,pass :: load => load_blocks_matrix
procedure,pass :: append => append_blocks_matrix
procedure,pass :: push => push_blocks_matrix
procedure,pass :: block => get_block_blocks_matrix
procedure,pass :: qn => get_qn_blocks_matrix
procedure,pass :: map => get_map_blocks_matrix
procedure,pass :: index => get_index_blocks_matrix
!
procedure,pass :: eigh => eigh_blocks_matrix
procedure,pass :: evals => evals_blocks_matrix
procedure,pass :: eval => eval_blocks_matrix
procedure,pass :: evec => evec_blocks_matrix
!
procedure,pass :: has_qn => has_qn_blocks_matrix
procedure,pass :: dump => dump_blocks_matrix
procedure,pass :: dgr => dgr_blocks_matrix
procedure,pass :: show => show_blocks_matrix
procedure,pass :: shape => shape_block_blocks_matrix
procedure,pass :: size => size_block_blocks_matrix
procedure,pass :: dims => dimensions_blocks_matrix
procedure,pass :: find => find_indices_blocks_matrix
procedure,pass :: sparse => sparse_blocks_matrix
end type blocks_matrix
interface blocks_matrix
module procedure :: construct_blocks_matrix
end interface blocks_matrix
interface as_blocks
module procedure :: construct_blocks_matrix
end interface as_blocks
!EQUALITY with scalar and function (A=B)
interface assignment(=)
module procedure :: blocks_matrix_equal_blocks_matrix
end interface assignment(=)
!INTRINSIC FUNCTION SIZE
intrinsic :: size
interface size
module procedure :: size_blocks_matrix
end interface size
!RETURN SHAPE OF THE WHOLE BLOCK MATRIX [Nrow,Ncol]
intrinsic :: shape
interface shape
module procedure :: shape_blocks_matrix
end interface shape
intrinsic :: transpose
interface transpose
module procedure :: transpose_blocks_matrix
end interface transpose
interface hconjg
module procedure :: hconjg_blocks_matrix
end interface hconjg
interface as_sparse
module procedure :: as_sparse_blocks_matrix
end interface as_sparse
! interface as_matrix
! module procedure :: as_matrix_blocks_matrix
! end interface as_matrix
public :: blocks_matrix
public :: as_blocks
public :: size
public :: shape
public :: transpose
public :: hconjg
public :: as_sparse
public :: assignment(=)
contains
!##################################################################
!##################################################################
! LIST CONSTRUCTOR/DESTRUCTOR
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: Intrinsic constructor
!+------------------------------------------------------------------+
function construct_blocks_matrix(matrix,qn,map) result(self)
real(8),dimension(:,:),intent(in) :: matrix
real(8),dimension(:),intent(in) :: qn
integer,dimension(:) :: map
type(blocks_matrix) :: self
call self%free()
allocate(self%root)
call self%append(matrix,qn,map)
end function construct_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: load a block matrix into a block one
!+------------------------------------------------------------------+
subroutine load_blocks_matrix(self,matrix,qn,map)
class(blocks_matrix),intent(inout) :: self
real(8),dimension(:) :: qn
integer,dimension(:) :: map
real(8),dimension(:,:),intent(in) :: matrix
call self%free()
allocate(self%root)
call self%append(matrix,qn,map)
end subroutine load_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: free an entire block matrix
!+------------------------------------------------------------------+
subroutine free_blocks_matrix(self)
class(blocks_matrix),intent(inout) :: self
type(block_type),pointer :: p,c
if(.not.associated(self%root))return
do
p => self%root
c => p%next
if(.not.associated(c))exit !empty list
p%next => c%next
c%next => null()
c%index= 0
if(allocated(c%qn))deallocate(c%qn)
if(allocated(c%M))deallocate(c%M)
if(allocated(c%map))deallocate(c%map)
deallocate(c)
enddo
if(allocated(self%evalues))deallocate(self%evalues)
if(allocated(self%eorder))deallocate(self%eorder)
self%Nblock=0
self%Nrow=0
self%Ncol=0
self%diag=.false.
self%root=>null()
p=>null()
c=>null()
end subroutine free_blocks_matrix
!##################################################################
!##################################################################
! APPEND/PUSH a block to the list
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: Append a block
!+------------------------------------------------------------------+
subroutine append_blocks_matrix(self,matrix,qn,map)
class(blocks_matrix),intent(inout) :: self
real(8),dimension(:,:),intent(in) :: matrix
real(8),dimension(:),intent(in) :: qn
integer,dimension(:) :: map
integer :: Dim
type(block_type),pointer :: p,c
!
if(.not.associated(self%root))allocate(self%root)
!
! Dim = size(matrix,1);call assert_shape(matrix,[Dim,Dim],"append_block","matrix")
!
p => self%root
c => p%next
do
if(.not.associated(c))exit
p => c
c => c%next
end do
!
allocate(p%next)
allocate(p%next%M, source=matrix)
allocate(p%next%map, source=map)
allocate(p%next%qn, source=qn)
p%next%index = p%index+1
!
if(.not.associated(c))then !end of the list special case (c=>c%next)
p%next%next => null()
else
p%next%next => c !the %next of the new node come to current
end if
self%Nblock = self%Nblock+1
!
self%Nrow = self%Nrow+size(Matrix,1)
self%Ncol = self%Ncol+size(Matrix,2)
!
p=>null()
c=>null()
end subroutine append_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: PUSH a block
!+------------------------------------------------------------------+
subroutine push_blocks_matrix(self,matrix,qn,map)
class(blocks_matrix),intent(inout) :: self
real(8),dimension(:,:),intent(in) :: matrix
real(8),dimension(:),intent(in) :: qn
integer,dimension(:) :: map
logical :: iupdate
type(block_type),pointer :: p,c
integer :: Dim
!
if(.not.associated(self%root))allocate(self%root)
!
Dim = size(matrix,1);call assert_shape(matrix,[Dim,Dim],"append_block","matrix")
!
if(.not.self%has_qn(qn))then
call self%append(matrix,qn,map)
return
endif
!
p => self%root
c => p%next
do
if(.not.associated(c))exit
if(all(c%qn == qn))then
exit
endif
p => c
c => c%next
end do
!
self%Nrow = self%Nrow-size(c%M,1)
self%Ncol = self%Ncol-size(c%M,2)
if(allocated(c%M))deallocate(c%M)
allocate(c%M, source=Matrix)
if(allocated(c%map))deallocate(c%map)
allocate(c%map, source=map)
if(allocated(c%qn))deallocate(c%qn)
allocate(c%qn, source=qn)
self%Nrow = self%Nrow+size(Matrix,1)
self%Ncol = self%Ncol+size(Matrix,2)
!
p=>null()
c=>null()
end subroutine push_blocks_matrix
!##################################################################
!##################################################################
! RETRIEVE CONTENT: GET BLOCK, DUMP WHOLE MATRIX
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: get
!+------------------------------------------------------------------+
function get_block_blocks_matrix(self,index,m) result(matrix)
class(blocks_matrix) :: self
integer,optional :: index,m
real(8),dimension(:,:),allocatable :: matrix
integer :: index_,m_
logical :: ifound
type(block_type),pointer :: c
!
index_=1
if(present(index))then
index_=index
elseif(present(m))then
m_ = m
if(self%diag)m_=self%eorder(m)
call self%find(m_,index_)
else
stop "get_block_blocks_matrix error: !present(index) + !present(m)"
endif
if(index_>self%Nblock.OR.index_<=0)stop "get_block_blocks_matrix error: block_index !in [1,self.size]"
!
ifound=.false.
c => self%root%next
do !traverse the list until QN is found
if(.not.associated(c))exit
if(c%index == index_) then
ifound=.true.
exit
endif
c => c%next
end do
if(.not.ifound)stop "get_blocks_matrix error: not found"
!
if(allocated(matrix))deallocate(matrix)
allocate(matrix, source=c%M)
! matrix = c%M
!
c=>null()
end function get_block_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: get
!+------------------------------------------------------------------+
function get_qn_blocks_matrix(self,index,m) result(qn)
class(blocks_matrix) :: self
integer,optional :: index
integer,optional :: m
integer :: index_,q,m_
real(8),dimension(:),allocatable :: qn
logical :: ifound
type(block_type),pointer :: c
!
if(allocated(qn))deallocate(qn)
!
index_=1
if(present(index))then
index_=index
elseif(present(m))then
m_ = m
if(self%diag)m_=self%eorder(m)
call self%find(m_,index_)
else
stop "get_qn_blocks_matrix error: !present(index) + !present(m)"
endif
if(index_>self%Nblock.OR.index_<=0)stop "get_qn_blocks_matrix error: block_index !in [1,self.size]"
ifound=.false.
c => self%root%next
do !traverse the list until QN is found
if(.not.associated(c))exit
if(c%index == index_) then
ifound=.true.
exit
endif
c => c%next
end do
if(.not.ifound)stop "get_qn_matrix error: not found"
!
allocate(qn, source=c%qn)
! qn = c%qn
!
c=>null()
end function get_qn_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: get
!+------------------------------------------------------------------+
function get_map_blocks_matrix(self,index,m) result(map)
class(blocks_matrix) :: self
integer,optional :: index
integer,optional :: m
integer :: index_,q,m_
integer,dimension(:),allocatable :: map
logical :: ifound
type(block_type),pointer :: c
!
index_=1
if(present(index))then
index_=index
elseif(present(m))then
m_ = m
if(self%diag)m_=self%eorder(m)
call self%find(m_,index_)
else
stop "get_qn_blocks_matrix error: !present(index) + !present(m)"
endif
if(index_>self%Nblock.OR.index_<=0)stop "get_qn_blocks_matrix error: block_index !in [1,self.size]"
ifound=.false.
c => self%root%next
do !traverse the list until QN is found
if(.not.associated(c))exit
if(c%index == index_) then
ifound=.true.
exit
endif
c => c%next
end do
if(.not.ifound)stop "get_qn_matrix error: not found"
!
if(allocated(map))deallocate(map)
allocate(map, mold=c%map)
map = c%map
!
c=>null()
end function get_map_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: get the element value at position (i,j) in the block matrix
!+------------------------------------------------------------------+
function get_index_blocks_matrix(self,qn) result(index)
class(blocks_matrix) :: self
real(8),dimension(:),intent(in) :: qn
integer :: index
logical :: ifound
type(block_type),pointer :: c
!
!qn does not exist: return with index=0
index=0
if(.not.self%has_qn(qn))return
!
!else qn exists so we can find it
c => self%root%next
do !traverse the list until QN is found
if(.not.associated(c))exit
if(all(c%qn == qn)) then
exit
endif
c => c%next
end do
!
index = c%index
!
c=>null()
end function get_index_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: dump a block matrix into a regular 2dim array
!+------------------------------------------------------------------+
function dump_blocks_matrix(self) result(matrix)
class(blocks_matrix),intent(inout) :: self
real(8),dimension(:,:),allocatable :: matrix,mtmp
integer :: Offset1,Offset2
integer :: N1,N2
type(block_type),pointer :: c
if(allocated(matrix))deallocate(matrix)
allocate(matrix(self%Nrow,self%Ncol));matrix=0d0
Offset1=0
Offset2=0
c => self%root%next
do
if(.not.associated(c))exit
N1 = size(c%M,1)
N2 = size(c%M,2)
Matrix(Offset1+1:Offset1+N1,Offset2+1:Offset2+N2) = c%M
Offset1 = Offset1 + N1
Offset2 = Offset2 + N2
c => c%next
end do
c=>null()
!
end function dump_blocks_matrix
! !+------------------------------------------------------------------+
! !PURPOSE: dump a sparse matrix into a regular 2dim array
! !+------------------------------------------------------------------+
! function as_matrix_blocks_matrix(self) result(matrix)
! class(blocks_matrix),intent(inout) :: self
! real(8),dimension(:,:),allocatable :: matrix
! if(allocated(matrix))deallocate(matrix)
! matrix = self%dump()
! end function as_matrix_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: dump a sparse matrix into a regular 2dim array
!+------------------------------------------------------------------+
function as_sparse_blocks_matrix(self) result(sparse)
class(blocks_matrix),intent(inout) :: self
type(sparse_matrix) :: sparse
integer,dimension(2) :: dims
integer :: i,it
real(8),dimension(:),allocatable :: self_vec
integer,dimension(:),allocatable :: self_map
dims = self%shape()
call sparse%init(dims(1),dims(2))
do it=1,dims(2)
self_vec = self%evec(m=it)
self_map = self%map(m=it)
do i=1,size(self_vec)
call sparse%insert(self_vec(i),self_map(i),it)
enddo
enddo
end function as_sparse_blocks_matrix
! !+------------------------------------------------------------------+
! !PURPOSE: dump a sparse matrix into a regular 2dim array
! !+------------------------------------------------------------------+
! function as_matrix_blocks_matrix(self) result(matrix)
! class(blocks_matrix),intent(inout) :: self
! type(sparse_matrix) :: sparse
! real(8),dimension(:,:),allocatable :: matrix
! if(allocated(matrix))deallocate(matrix)
! matrix = as_matrix(as_sparse(self))
! end function as_matrix_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: dump a sparse matrix into a regular 2dim array
!+------------------------------------------------------------------+
function sparse_blocks_matrix(self,n,m) result(sparse)
class(blocks_matrix),intent(inout) :: self
integer :: n,m
type(sparse_matrix) :: sparse
integer,dimension(2) :: dims
integer :: i,it
real(8),dimension(:),allocatable :: self_vec
integer,dimension(:),allocatable :: self_map
! dims = shape(self)
! m_=dims(2);if(present(m))m_=m
! if(m_<1.OR.m_>dims(2))stop "as_sparse_truncate_blocks_matrix ERROR: m<1 OR m>size(self,2)"
call sparse%init(n,m)
do it=1,m
self_vec = self%evec(m=it)
self_map = self%map(m=it)
do i=1,size(self_vec)
call sparse%fast_insert(self_vec(i),self_map(i),it)
enddo
enddo
end function sparse_blocks_matrix
!##################################################################
!##################################################################
! LINEAR ALGEBRA
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: performs eigh on each block of the blocks_matrix,
! the EigenVectors matrix overwrites the block and the EigenValues
! are stored in the array E of the block (so far unused).
!+------------------------------------------------------------------+
subroutine eigh_blocks_matrix(self,sort,reverse,order)
class(blocks_matrix),intent(inout) :: self
logical,optional :: sort,reverse
integer,dimension(:),allocatable,optional :: order
real(8),dimension(:),allocatable :: Rtmp
integer,dimension(:),allocatable :: Itmp
logical :: sort_,reverse_
type(block_type),pointer :: c
integer :: i,Nloc,Offset,N
!
sort_ =.true.;if(present(sort))sort_=sort
reverse_=.true.;if(present(reverse))reverse_=reverse
!
N = self%Nrow
if(N/=self%Ncol)print*,"WARNING eigh blocks matrix: self is not square"
if(allocated(self%evalues))deallocate(self%evalues)
if(allocated(self%eorder))deallocate(self%eorder)
allocate(self%evalues(N));self%evalues=0d0
allocate(self%eorder(N));self%eorder=0
!
Offset=0
c => self%root%next
do !loop over all blocks
if(.not.associated(c))exit
Nloc = size(c%M,1)
if(any(shape(c%M)/=[Nloc,Nloc]))stop "eigh block matrix ERROR: local block is not square"
if(allocated(c%E))deallocate(c%E)
allocate(c%E(Nloc));c%E=0d0
call eigh(c%M,c%E) !<- overwrites blocks with eigenvec matrix
!
self%evalues(Offset+1:Offset+Nloc) = c%E
Offset = Offset + Nloc
!
c => c%next
end do
c=>null()
!
!init Eorder vector to 1..N
self%eorder=(/(i,i=1,N)/)
if(sort_)then
call sort_quicksort(self%evalues,self%eorder)
if(reverse_)then
allocate(Rtmp(N),Itmp(N))
Rtmp = self%evalues(N:1:-1) ; self%evalues=Rtmp
Itmp = self%eorder(N:1:-1) ; self%eorder =Itmp
deallocate(Rtmp,Itmp) !just to avoid in place reverse
endif
endif
if(present(order))then
if(allocated(order))deallocate(order)
allocate(order, source=self%eorder)
endif
!
self%diag=.true.
!
end subroutine eigh_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: returns an array with all evals of the block matrix
!+------------------------------------------------------------------+
function evals_blocks_matrix(self) result(evals)
class(blocks_matrix) :: self
real(8),dimension(:),allocatable :: evals
type(block_type),pointer :: c
!
if(.not.self%diag)stop "Evals_blocks_matrix error: self.diag=F, call self.eigh() before trying to get an eigenvector"
!
if(allocated(evals))deallocate(evals)
!
allocate(evals, source=self%evalues)
!
end function evals_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: returns a single eigenvalue with all evals of the block matrix
!+------------------------------------------------------------------+
function eval_blocks_matrix(self,m) result(eval)
class(blocks_matrix) :: self
integer :: m
real(8) :: eval
if(.not.self%diag)stop "Evals_blocks_matrix error: self.diag=F, call self.eigh() before trying to get an eigenvector"
if(m>self%Nrow.OR.m<=0)stop "eval_block_matrix warning: m !in [1,self.Ndim]"
eval = self%evalues(m)
end function eval_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE:
!+------------------------------------------------------------------+
function evec_blocks_matrix(self,m) result(vec)
class(blocks_matrix) :: self
integer :: m
real(8),dimension(:),allocatable :: vec
integer :: m_,i,q,pos
type(block_type),pointer :: c
!
! if(.not.self%diag)stop "Evec_blocks_matrix error: self.diag=F, call self.eigh() before trying to get an eigenvector"
!
if(m>self%Ncol.OR.m<=0)stop "evec_block_matrix warning: m !in [1,self.Ncol]"
m_=m
if(self%diag)m_=self%eorder(m)
call self%find(m_,q,pos)
!
c => self%root
do i=1,q
c => c%next
enddo
!
if(allocated(vec))deallocate(vec)
allocate(vec, source=c%M(:,pos))
!
c=>null()
end function evec_blocks_matrix
!##################################################################
!##################################################################
! OVERLOAD INTRINSIC and OPERATIONS
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: Returns the size of given qmap
!+------------------------------------------------------------------+
function size_blocks_matrix(self) result(size)
class(blocks_matrix),intent(in) :: self
integer :: size
size = self%Nblock
end function size_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: Return the shape of a sparse matrix
!+------------------------------------------------------------------+
function shape_blocks_matrix(self) result(shape)
class(blocks_matrix),intent(in) :: self
integer,dimension(2) :: shape
shape = [self%Nrow,self%Ncol]
end function shape_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE:
!+------------------------------------------------------------------+
function shape_block_blocks_matrix(self,index) result(Bshape)
class(blocks_matrix) :: self
integer,optional :: index
integer :: index_
integer,dimension(2) :: Bshape
logical :: ifound
type(block_type),pointer :: c
!
if(.not.present(index))then
Bshape = [self%Nrow,self%Ncol]
return
endif
!
index_=1;if(present(index))index_=index
if(index_>self%Nblock.OR.index_<=0)stop "shape_block_blocks_matrix error: block_index !in [1,self.size]"
!
ifound=.false.
c => self%root%next
do !traverse the list until QN is found
if(.not.associated(c))exit
if(c%index == index_) then
ifound=.true.
exit
endif
c => c%next
end do
if(.not.ifound)stop "get_blocks_matrix error: not found"
!
Bshape = shape(c%M)
!
c=>null()
end function shape_block_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE:
!+------------------------------------------------------------------+
function size_block_blocks_matrix(self,index,dim) result(Bsize)
class(blocks_matrix) :: self
integer,optional :: index
integer,optional :: dim
integer :: index_
integer :: dim_
integer :: Bsize
integer :: Mshape(2)
logical :: ifound
type(block_type),pointer :: c
!
dim_=1;if(present(dim))dim_=dim
if(dim_ < 1 .OR. dim_ > 2)then
write(0,*)"get warning: dim_ != [1,2]: reset to dim=1"
dim_=1
endif
!
index_=1;if(present(index))index_=index
if(index_>self%Nblock.OR.index_<=0)stop "shape_block_blocks_matrix error: block_index !in [1,self.size]"
!
if(.not.present(index))then !return the dimension of the whole self.
Mshape= [self%Nrow,self%Ncol]
Bsize = Mshape(dim_)
return
endif
!
ifound=.false.
c => self%root%next
do !traverse the list until QN is found
if(.not.associated(c))exit
if(c%index == index_) then
ifound=.true.
exit
endif
c => c%next
end do
if(.not.ifound)stop "get_blocks_matrix error: not found"
!
Bsize = size(c%M,dim_)
!
c=>null()
end function size_block_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: Returns dimensions D_q=1,size(self) of the Block matrix
!+------------------------------------------------------------------+
function dimensions_blocks_matrix(self,dim) result(Dvec)
class(blocks_matrix),intent(inout) :: self
integer,dimension(:),allocatable :: Dvec
integer,optional :: dim
integer :: dim_
type(block_type),pointer :: c
integer :: q
dim_=1;if(present(dim))dim_=dim
if(dim_ < 1 .OR. dim_ > 2)then
write(0,*)"get warning: dim_ != [1,2]: reset to dim=1"
dim_=1
endif
if(allocated(Dvec))deallocate(Dvec)
allocate(Dvec(size(self)))
do q=1,size(self)
Dvec(q) = self%size(index=q,dim=dim_)
enddo
end function dimensions_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: Given an integer index m=1,self.Nrow,
! returns the corresponding BLock index containing it and the relative
! position pos inside the block.
! [1...D_1][D_1+1...D_2]...[D_{q-1}...D_q] = D=self.Ndim
!+------------------------------------------------------------------+
subroutine find_indices_blocks_matrix(self,m,index,pos)
class(blocks_matrix) :: self
integer :: m
integer :: index
integer :: q
integer,optional :: pos
integer :: pos_
integer,dimension(:),allocatable :: Dq
Dq = self%dims(dim=2)
if(m>sum(Dq).OR.m<=0)stop "find_indices_blocks_matrix error: m !in [1,D==self.Ncol]"
do q=1,size(Dq)
if(m <= sum(Dq(1:q)))then
index = q
pos_ = m - sum(Dq(1:q-1))
exit
endif
enddo
if(present(pos))pos=pos_
end subroutine find_indices_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: Returns True is qn exists, False otherwise
!+------------------------------------------------------------------+
function has_qn_blocks_matrix(self, qn) result(bool)
class(blocks_matrix),intent(inout) :: self
real(8),dimension(:),intent(in) :: qn
type(block_type),pointer :: c
logical :: bool
bool=.false.
c => self%root%next
do !traverse the list until index is found
if(.not.associated(c))exit
if(all(c%qn == qn)) then
bool=.true.
exit
endif
c => c%next
end do
end function has_qn_blocks_matrix
!##################################################################
!##################################################################
! SHOW
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: show matrix
!+------------------------------------------------------------------+
subroutine show_blocks_matrix(self,fmt,file)
class(blocks_matrix),intent(in) :: self
character(len=*),optional :: fmt
character(len=*),optional :: file
character(len=12) :: fmt_
integer :: i,j,unit_
character(len=64) :: format
real(8) :: val
type(block_type),pointer :: c
!
unit_=6
fmt_=str(show_fmt);if(present(fmt))fmt_=str(fmt)
if(present(file))open(free_unit(unit_),file=str(file))
!
format='('//str(fmt_)//'1x)'
! format='(A1,'//str(fmt_)//',A1,'//str(fmt_)//',A1,1x)'
!
write(unit_,"(A6,I12)")"Size :",size(self)
write(unit_,"(A6,2I6)")"Shape:",shape(self)
write(unit_,"(A18)")"------------------"
c => self%root%next
do
if(.not.associated(c))exit
write(unit_,"(A6,I12)")"Index:",c%index
write(unit_,"(A6,"//str(size(c%map))//"I6)")"map :",c%map
write(unit_,"(A6,"//str(size(c%qn))//"F12.4)")"QN :",c%qn
if(allocated(c%E))&
write(unit_,"(A6,"//str(size(c%E))//"F12.4)")"E :",c%E
write(unit_,"(A6)")"Block:"
do i=1,size(c%M,1)
do j=1,size(c%M,2)
val = c%M(i,j)
write(unit_,"("//str(self%Ncol)//"(F12.4,1X))",advance='no')val
!write(unit_,"("//str(self%Ncol)//str(format)//")",advance='no')"(",dreal(val),",",dimag(val),")"
enddo
write(unit_,*)
enddo
write(unit_,*)
c => c%next
enddo
if(present(file))close(unit_)
end subroutine show_blocks_matrix
!##################################################################
!##################################################################
! BLOCK MATRIX BASIC ALGEBRA
!##################################################################
!##################################################################
function dgr_blocks_matrix(a) result(adg)
class(blocks_matrix), intent(in) :: a
type(blocks_matrix) :: adg
integer :: i
call adg%free()
do i=1,size(a)
call adg%append( (transpose(a%block(index=i))), a%qn(index=i), a%map(index=i))
enddo
end function dgr_blocks_matrix
function transpose_blocks_matrix(a) result(adg)
class(blocks_matrix), intent(in) :: a
type(blocks_matrix) :: adg
integer :: i
call adg%free()
do i=1,size(a)
call adg%append((transpose(a%block(index=i))), a%qn(index=i), a%map(index=i))
enddo
end function transpose_blocks_matrix
function hconjg_blocks_matrix(a) result(adg)
class(blocks_matrix), intent(in) :: a
type(blocks_matrix) :: adg
integer :: i
call adg%free()
do i=1,size(a)
call adg%append((transpose(a%block(index=i))), a%qn(index=i), a%map(index=i))
enddo
end function hconjg_blocks_matrix
!+------------------------------------------------------------------+
!PURPOSE: Block matrix equality A = B. Deep copy
!+------------------------------------------------------------------+
subroutine blocks_matrix_equal_blocks_matrix(a,b)
type(blocks_matrix),intent(inout) :: a
type(blocks_matrix),intent(in) :: b
integer :: i
call a%free()
do i=1,size(b)
call a%append(b%block(index=i), b%qn(index=i), b%map(index=i))
enddo
end subroutine blocks_matrix_equal_blocks_matrix
end module MATRIX_BLOCKS
!##################################################################
!##################################################################
!##################################################################
!##################################################################
! /_ __/ ____/ ___/_ __/
! / / / __/ \__ \ / /
! / / / /___ ___/ // /
! /_/ /_____//____//_/
!##################################################################
!##################################################################
!##################################################################
!##################################################################
#ifdef _TEST
program testBLOCK_MATRICES
USE SCIFOR
USE MATRIX_BLOCKS
USE TUPLE_BASIS
USE LIST_SECTORS
implicit none
type(blocks_matrix) :: blH,a,adg,bvec(3)
type(tbasis) :: a_basis
type(sectors_list) :: a_sector