-
Notifications
You must be signed in to change notification settings - Fork 20
/
group_mod.f90
1775 lines (1544 loc) · 54.8 KB
/
group_mod.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
!-----------------------------------------------------------------------
! PSCF - Polymer Self-Consistent Field Theory
! Copyright (2002-2016) Regents of the University of Minnesota
! contact: David Morse, [email protected]
!
! This program is free software; you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation. A copy of this license is included in
! the LICENSE file in the top-level PSCF directory.
!-----------------------------------------------------------------------
!****m scf/group_mod
! MODULE
! group_mod
! PURPOSE
! Define derived data types and basic operations
! for space group symmetry operations and groups
! AUTHOR
! David Morse (2002)
! SOURCE
!-----------------------------------------------------------------------
module group_mod
use const_mod, only : dim, long
use version_mod, only : version_type, input_version, output_version
implicit none
private
! Derived types
public :: symmetry_type ! space group symmetry
public :: group_type ! space group
! Generic interfaces
public :: operator(.dot.) ! .dot. products for vectors,
! matrices and symmetry_type
public :: inverse ! inversion of 2D and 3D real
! matrices and symmetry_type
public :: equal ! equality with tolerance for variety
! of data types, including symmetries
! Public procedures
public :: read_group, output_group ! io for groups
public :: make_group ! complete and check space group
!***
!Private (make public for testing)
!public :: max_order
!public :: output_symmetry
!public :: table_type, make_table
!public :: valid_basis
!public :: index_of_inverse
!public :: make_identity, shift_translation
!public :: is_sub_group
! Parameters (Private)
integer, parameter :: max_order = 192
real(long), parameter :: epsilon = 1.0E-6_long
character(9), parameter :: Cartesian = 'Cartesian'
character(9), parameter :: Bravais = 'Bravais '
!------------------------------------------------------------------
!****t group_mod/symmetry_type
! TYPE
! symmetry_type
! VARIABLE
! character(9) basis = 'Cartesian' or 'Bravais '
! real(long) m(3,3) = point group matrix
! real(long) v(3) = translation vector
! COMMENT
! Conventions for symmetry_type and related derived types:
!
! a) The effect of a symmetry on a position vector is to take
!
! R -> m .dot. R + v
!
! where m .dot. R represents contraction with first index of m
!
! b) Point group matrix m operates on reciprocal G vectors by
! contraction with first index of m: G -> G .dot. m
!
! c) Symmetries can be expressed in either Cartesian or
! Bravais basis, as indicated by value of the character
! variable basis, which can have values equal to the string
! constants Cartesian = 'Cartesian' or Bravais='Bravais '.
!
! d) In the Bravais basis, position vectors are represented
! as coefficients in expansion in Bravais basis vectors,
! R_basis(:,i), i=1,..dim, and G vectors are represented
! as (integer) coefficients in expansion in reciprocal
! basis vectors, G_basis(:,j), j=1,..,dim .
!
! e) In the Bravais representation, elements of a point
! group matrix m should be integers (though they are
! stored as reals), and elements of the translation
! vector v should be low order fractions
!
! SOURCE
!------------------------------------------------------------------
type symmetry_type
character(9) :: basis ! must equal 'Cartesian' or 'Bravais '
real(long) :: m(3,3) ! point group matrix
real(long) :: v(3) ! translation vector
end type symmetry_type
!***
!------------------------------------------------------------------
!****t group_mod/group_type
! TYPE
! group_type
! VARIABLE
! order = # of symmetry elements in group
! s(max_order) = array of symmetries
! COMMENT
! All symmetries in a group must have the same basis, i.e.,
! they must all be either in Bravais or Cartesian basis
! SOURCE
!------------------------------------------------------------------
type group_type
integer :: order
type(symmetry_type) :: s(max_order)
end type group_type
!***
!------------------------------------------------------------------
!****it group_mod/table_type
! TYPE
! table_type - multipication table for symmetries in a group
!
! VARIABLE
!
! order = order of symmetries
! s(i,j) = product symmetries i .dot. symmetry j
! index(i,j) = index of s(i,j) in group
! so that s(i,j) = group%s( index(i,j) )
!
! SOURCE
!------------------------------------------------------------------
type table_type
integer :: order
integer :: index(max_order,max_order)
type(symmetry_type) :: s(max_order,max_order)
end type table_type
!***
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
! Generic Interfaces
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!------------------------------------------------------------------
!****p group_mod/dot
! OPERATOR
! .dot. - dot product generic interface
!
! COMMENT
!
! a) In names of the specific realizations the types of
! arguments are indicated with the shorthand:
!
! integer ivec(:)
! real(long) vec(:)
! real(long) mat(:,:)
! symmetry_type sym
!
!
! b) When evaluating dot products, elements of the input
! arguments with indices > dim are ignored, and elements
! of any returned vector or matrix with indices > dim
! are padded with zeros. Because the return value of
! the operator cannot be adjustable, the operator returns
! vectors and matrices (when appropriate) with dimension 3
!
! SOURCE
!------------------------------------------------------------------
interface operator(.dot.)
module procedure ivec_dot_ivec ! integer
module procedure vec_dot_vec ! real
module procedure ivec_dot_vec ! real
module procedure vec_dot_ivec ! real
module procedure mat_dot_vec ! real(3)
module procedure ivec_dot_mat ! real(3)
module procedure vec_dot_mat ! real(3)
module procedure mat_dot_mat ! real(3,3)
module procedure sym_dot_vec ! real(3) = sym%m.dot.vec + sym%v
module procedure vec_dot_sym ! real(3) = vec.dot.sym%m
module procedure ivec_dot_sym ! integer(3) = ivec.dot.sym%m
module procedure sym_dot_sym ! symmetry_type
end interface
!***
!------------------------------------------------------------------
!****p group_mod/equal
! FUNCTION
! equal(a,b)
! RETURN
!
! true if a nd b are equal to within a tolerance epsilon
!
! The arguments a and b may be objects of type:
!
! real(long)
! integer dimension(3)
! real(long) dimension(3)
! real(long) dimension(3,3)
! symmetry_type
!
! SOURCE
!------------------------------------------------------------------
interface equal
module procedure real_equal
module procedure ivector_equal
module procedure r_vector_equal
module procedure matrix_equal
module procedure symmetry_equal
end interface
!***
!------------------------------------------------------------------
!****p group_mod/inverse
! FUNCTION
!
! inverse(a) - generic interface for inversion
!
! RETURN
!
! inverse of argument a
!
! The argument a may be:
!
! real(long) matrix_inverse(3,3)
! symmetry_type symmetry_inverse
!
! SOURCE
!------------------------------------------------------------------
interface inverse
module procedure matrix_inverse ! (3,3) padded with zeros
module procedure symmetry_inverse ! symmetry_type
end interface
!***
contains
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
! Output Subroutines for Data Types & Input for Groups
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!------------------------------------------------------------------
!****ip group_mod/output_scalar
! SUBROUTINE
! output_scalar(r,iunit)
! PURPOSE
! Write real scalar to iunit
! SOURCE
!------------------------------------------------------------------
subroutine output_scalar(r,iunit)
real(long), intent(IN) :: r
integer :: iunit
!***
write(iunit,FMT='(F12.3)') r
write(iunit,*)
end subroutine output_scalar
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/output_vector
! SUBROUTINE
! output_vector(v,iunit)
! PURPOSE
! Write real vector to iunit
! SOURCE
!------------------------------------------------------------------
subroutine output_vector(v,iunit)
real(long), intent(in) :: v(3)
integer, intent(in) :: iunit
!***
integer :: i
! Check vector size
if (size(v) < dim) then
write(6,*) 'Error: Too small a vector in output_matrix'
endif
write(iunit,FMT='(3F12.3)') (v(i), i=1, dim)
write(iunit,*)
end subroutine output_vector
!===================================================================
!-------------------------------------------------------------------
!****ip group_mod/output_matrix
! SUBROUTINE
! output_matrix(m,iunit)
! PURPOSE
! Write matrix m to file iunit
! SOURCE
!-------------------------------------------------------------------
subroutine output_matrix(m,iunit)
real(long), intent(in) :: m(:,:)
integer :: iunit
!***
integer :: i, j
! Check matrix size
if ( (size(m,1) < dim).or.(size(m,2)< dim)) then
write(6,*) 'Error: Too small a matrix in output_matrix'
endif
do i=1, dim
write(iunit,FMT='(3F12.3)') (m(i,j), j=1, dim)
enddo
write(iunit,*)
end subroutine output_matrix
!===================================================================
!-------------------------------------------------------------------
!****p group_mod/output_symmetry
! SUBROUTINE
! output_symmetry(s, iunit, mode_arg)
! PURPOSE
! Write symmetry s to file iunit
! SOURCE
!-------------------------------------------------------------------
subroutine output_symmetry(s, iunit, mode_arg)
type(symmetry_type), intent(in) :: s
integer :: iunit
integer, intent(in), optional :: mode_arg
!***
integer :: i, j, mode
! Set output format mode
if (present(mode_arg)) then
mode = mode_arg
else
mode = 1
end if
if (mode == 1) then
! Output in native (fortran pscf) format
do i=1, dim
if (dim == 3) then
write(iunit,FMT='(3F8.3,5X,F8.3)') &
(s%m(i,j),j=1, dim), s%v(i)
else if (dim == 2) then
write(iunit,FMT='(2F8.3,5X,F8.3)') &
(s%m(i,j),j=1, dim), s%v(i)
else if (dim == 1) then
write(iunit,FMT='(F8.3,5X,F8.3)') &
(s%m(i,j),j=1, dim), s%v(i)
else
write(iunit,*) 'Invalid dim in output_symmetry'
endif
enddo
write(iunit,*)
else if (mode == 2) then
! Output in pscfpp (c++ code) format
call output_symmetry_c(s, iunit);
end if
end subroutine output_symmetry
!===================================================================
!-------------------------------------------------------------------
!****p group_mod/output_symmetry_c
! SUBROUTINE
! output_symmetry_c(s,iunit)
! PURPOSE
! Write symmetry s to file iunit in format of C++ pscfpp code.
! SOURCE
!-------------------------------------------------------------------
subroutine output_symmetry_c(s, iunit)
type(symmetry_type), intent(in) :: s
integer :: iunit
!***
real(long) :: r, q
integer :: i, j, k, d, n
integer :: rot(3,3) ! point group matrix
integer :: den(3) ! denominator of translation fraction
integer :: num(3) ! numerator of translation fraction
integer :: trial(4) ! allowed values of denominator
logical :: found
! Convert elements of matrix m to nearest integers
do i=1, dim
do j=1, dim
r = anint(s%m(i, j))
if (dabs(r - s%m(i,j)) > epsilon) then
write(iunit,*) 'Rotation matrix element is non-integer'
return
else
rot(i, j) = nint(r)
endif
enddo
enddo
! Create allowed values of denominator
trial(1) = 2
trial(2) = 3
trial(3) = 4
trial(4) = 6
! Convert elements of vector s%v to nearest fractions
do i=1, dim
found = .false.
r = anint(s%v(i))
if (abs(r - s%v(i)) < 1.0E-3) then
found = .true.
num(i) = 0
den(i) = 1
else
do j=1, 4
k = trial(j)
q = k*s%v(i)
r = anint(q)
if (abs(r - q) < 1.0E-3) then
found = .true.
den(i) = k
num(i) = nint(q)
q = real(num(i))/real(den(i))
if (abs(q - s%v(i)) > 1.0E-2) then
write(iunit,*) 'Incorrect conversion to fraction'
write(iunit,*) 'Translation element = ', s%v(i)
write(iunit,*) 'Numerator = ', num(i)
write(iunit,*) 'Denominator = ', den(i)
return
endif
exit
endif
enddo
if (den(i) < 2) then
write(iunit,*) 'Invalid denominator'
write(iunit,*) 'Translation element = ', s%v(i)
write(iunit,*) 'Numerator = ', num(i)
write(iunit,*) 'Denominator = ', den(i)
return
endif
endif
if (.not.found) then
write(iunit,*) 'Translation element = ', s%v(i)
write(iunit,*) 'Translation element is not small fraction'
return
endif
enddo
! Write rotation matrix
do i=1, dim
if (dim == 3) then
write(iunit,FMT='(3I3)') (rot(i,j), j=1, dim)
else if (dim == 2) then
write(iunit,FMT='(2I3)') (rot(i,j), j=1, dim)
else if (dim == 1) then
write(iunit,FMT='(I3)') (rot(i,j), j=1, dim)
else
write(iunit,*) 'Invalid dim in output_symmetry'
endif
enddo
! Write translation vector
do i=1, dim
if (num(i) == 0) then
write(iunit,FMT='(I3)', advance="no") 0
else
write(iunit,FMT='(I3,"/",I1)', advance="no") num(i), den(i)
endif
enddo
write(iunit,*)
write(iunit,*)
end subroutine output_symmetry_c
!===================================================================
!-------------------------------------------------------------------
!****p group_mod/output_group
! SUBROUTINE
! output_group(g,iunit,mode_arg)
! PURPOSE
! Write group g to file iunit
! SOURCE
!-------------------------------------------------------------------
subroutine output_group(g,iunit,mode_arg)
type(group_type), intent(IN) :: g
integer, intent(IN) :: iunit
integer, intent(IN), optional :: mode_arg
!***
type(version_type) :: version
integer :: mode
integer :: i
! Check for optional mode argument
if (present(mode_arg)) then
mode = mode_arg
else
mode = 1
end if
! Output preamble (if any)
if (mode == 1) then
! Write preamble
version%major = 1
version%minor = 0
call output_version(version, iunit)
write(iunit,FMT="(A9,' = basis')") g%s(1)%basis
write(iunit,FMT="(i9,' = dimension')") dim
write(iunit,FMT="(i9,' = # of symmetries')") g%order
write(iunit,*)
! Write symmetry elements
do i=1, g%order
write(iunit,FMT='(i3)') i
call output_symmetry(g%s(i),iunit, mode)
enddo
else if (mode == 2) then
write(iunit,FMT="('dim ', i3)") dim
write(iunit,FMT="('size ', i3)") g%order
write(iunit,*)
! Write symmetry elements
do i=1, g%order
call output_symmetry(g%s(i),iunit, mode)
enddo
else
write(iunit, *) "Incorrect mode format number = ", mode
endif
end subroutine output_group
!===================================================================
!-------------------------------------------------------------------
!****ip group_mod/output_table
! SUBROUTINE
! output_table(g,iunit)
! PURPOSE
! Write table t to file iunit
! SOURCE
!-------------------------------------------------------------------
subroutine output_table(t,iunit)
type(table_type), intent(IN) :: t
integer, intent(IN) :: iunit
!***
integer :: i,j
do i=1, t%order
write(iunit,FMT='(48I3)') (t%index(i,j),j=1,t%order)
enddo
end subroutine output_table
!===================================================================
!-------------------------------------------------------------------
!****p group_mod/read_group
! SUBROUTINE
! read_group(g,iunit)
! PURPOSE
! Read group g from file unit iunit
! SOURCE
!-------------------------------------------------------------------
subroutine read_group(g,iunit)
type(group_type), intent(OUT) :: g
integer, intent(IN) :: iunit
!***
integer :: i, j, k, i_input, dim_input
character(9) :: basis
type(version_type) :: version
! Read file format
call input_version(version,iunit)
! Read character string basis and check validity
read(iunit,*) basis
if ( ( basis /= Cartesian ).and.( basis /= Bravais ) ) then
write(6,*) 'Invalid value of s%basis in read_group'
stop
endif
! Read dimension, and check consistency
read(iunit,*) dim_input
if ( dim_input /= dim ) then
write(6,*) 'Incompatible values of dim in read_group'
write(6,*) 'dim_input = ', dim_input
write(6,*) 'dim = ', dim
stop
endif
do i=1, max_order
g%s(i)%m = 0.0_long
g%s(i)%v = 0.0_long
enddo
read(iunit,*) g%order
do i=1, g%order
read(iunit,*)
read(iunit,*) i_input
do j=1, dim
read(iunit,*) (g%s(i)%m(j,k),k=1, dim),g%s(i)%v(j)
enddo
g%s(i)%basis = basis
enddo
end subroutine read_group
!===================================================================
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
! Definitions of Dot Product
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!------------------------------------------------------------------
!****ip group_mod/ivec_dot_ivec
! FUNCTION
! integer ivec_dot_ivec(v1,v2)
! RETURN
! Dot product of 2 integer vectors
! SOURCE
!------------------------------------------------------------------
integer function ivec_dot_ivec(v1,v2)
integer, intent(IN) :: v1(:), v2(:)
!***
integer :: i
ivec_dot_ivec = 0
do i=1, dim
ivec_dot_ivec = ivec_dot_ivec + v1(i)*v2(i)
enddo
end function ivec_dot_ivec
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/vec_dot_vec
! FUNCTION
! real(long) vec_dot_vec(v1,v2)
! RETURN
! Dot product of 2 real vectors
! SOURCE
!------------------------------------------------------------------
real(long) function vec_dot_vec(v1,v2)
real(long),intent(IN), dimension(:) :: v1, v2
!***
integer ::i
vec_dot_vec = 0.0_long
do i=1, dim
vec_dot_vec = vec_dot_vec + v1(i)*v2(i)
enddo
end function vec_dot_vec
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/ivec_dot_vec
! FUNCTION
! real(long) ivec_dot_vec(v1,v2)
! RETURN
! Dot product of integer vector .dot. real vector
! SOURCE
!------------------------------------------------------------------
real(long) function ivec_dot_vec(v1,v2)
integer, intent(IN), dimension(:) :: v1
real(long), intent(IN), dimension(:) :: v2
!***
integer ::i
ivec_dot_vec = 0.0_long
do i=1, dim
ivec_dot_vec = ivec_dot_vec + v1(i)*v2(i)
enddo
end function ivec_dot_vec
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/vec_dot_ivec
! FUNCTION
! real(long) vec_dot_ivec(v1,v2)
! RETURN
! Dot product of real vector .dot. integer vector
! SOURCE
!------------------------------------------------------------------
real(long) function vec_dot_ivec(v1,v2)
real(long), intent(IN), dimension(:) :: v1
integer, intent(IN), dimension(:) :: v2
!***
vec_dot_ivec = ivec_dot_vec(v2,v1)
end function vec_dot_ivec
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/mat_dot_vec
! FUNCTION
! real(long) mat_dot_vec(m,v)
! RETURN
! Dot product of real matrix m .dot. real column vector v
! Returns real array mat_dot_vec(3)
! SOURCE
!------------------------------------------------------------------
function mat_dot_vec(m,v)
real(long) :: mat_dot_vec(3)
real(long),intent(IN) :: v(:)
real(long),intent(IN) :: m(:,:)
!***
integer ::i,j
mat_dot_vec = 0.0_long
do i=1, dim
do j=1, dim
mat_dot_vec(i) = mat_dot_vec(i) + m(i,j)*v(j)
enddo
enddo
end function mat_dot_vec
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/vec_dot_mat
! FUNCTION
! real(long) vec_dot_mat(v,m)
! RETURN
! Dot product of real row vector v .dot. real matrix m
! Returns real array vec_dot_mat(3)
! SOURCE
!------------------------------------------------------------------
function vec_dot_mat(v,m)
real(long) :: vec_dot_mat(3)
real(long),intent(IN) :: v(:)
real(long),intent(IN) :: m(:,:)
!***
integer ::i,j
vec_dot_mat = 0.0_long
do j=1, dim
do i=1, dim
vec_dot_mat(j) = vec_dot_mat(j) + v(i)*m(i,j)
enddo
enddo
end function vec_dot_mat
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/ivec_dot_mat
! FUNCTION
! real(long) ivec_dot_mat(v,m)
! RETURN
! Dot product of integer row vector v .dot. real matrix m
! Returns real array ivec_dot_mat(3)
! SOURCE
!------------------------------------------------------------------
function ivec_dot_mat(v,m)
real(long) :: ivec_dot_mat(3)
integer, intent(IN) :: v(:)
real(long), intent(IN) :: m(:,:)
!***
integer ::i,j
ivec_dot_mat = 0.0_long
do j=1, dim
do i=1, dim
ivec_dot_mat(j) = ivec_dot_mat(j) + v(i)*m(i,j)
enddo
enddo
end function ivec_dot_mat
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/mat_dot_mat
! FUNCTION
! real(long) mat_dot_mat(m1,m2)
! RETURN
! Dot product of real matrix m1 .dot. real matrix m2
! Returns real array mat_dot_mat(3,3)
! SOURCE
!------------------------------------------------------------------
function mat_dot_mat(m1,m2)
real(long) :: mat_dot_mat(3,3)
real(long),intent(IN) :: m1(:,:), m2(:,:)
!***
integer ::i,j,k
mat_dot_mat = 0.0_long
do i=1, dim
do j=1, dim
do k=1, dim
mat_dot_mat(i,j) = mat_dot_mat(i,j) + m1(i,k)*m2(k,j)
enddo
enddo
enddo
end function mat_dot_mat
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/sym_dot_vec
! FUNCTION
! real(long) sym_dot_vec(s,v)
! RETURN
! Result of symmetry element s acting on position vector v
! Returns real vector sym_dot_vec(3) = s%m .dot. v + s%v
! COMMENT
! Result is valid iff basis of symmetry and vector are same
! (i.e., both Cartesian or both Bravais bases)
! SOURCE
!-------------------------------------------------------------------
function sym_dot_vec(s,v)
real(long), dimension(3) :: sym_dot_vec
type(symmetry_type), intent(IN):: s
real(long),intent(IN) :: v(:)
!***
integer:: i
do i=1, dim
sym_dot_vec = mat_dot_vec(s%m,v) + s%v
enddo
end function sym_dot_vec
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/vec_dot_sym
! FUNCTION
! real(long) vec_dot_sym(v,s)
! RETURN
! Result of symmetry s upon real Cartesian wavevector v
! Returns real array vec_dot_sym(3) = v .dot. s%m
! COMMENT
! Valid only if symmetry and wavector are in Cartesian basis
! SOURCE
!-------------------------------------------------------------------
function vec_dot_sym(v,s)
real(long), dimension(3) :: vec_dot_sym
real(long), intent(IN), dimension(:) :: v
type(symmetry_type), intent(IN) :: s
!***
if ( s%basis == Cartesian ) then
vec_dot_sym = vec_dot_mat(v,s%m)
else
write(6,*) 'Improper use of vec_dot_sym - not a Cartesian symmetry'
stop
endif
end function vec_dot_sym
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/ivec_dot_sym
! FUNCTION
! real(long) ivec_dot_sym(v,s)
! RETURN
! Result of symmetry s upon integer Cartesian wavevector v
! Returns integer array vec_dot_sym(3) = v .dot. s%m
! COMMENT
! Valid only if symmetry and wavector are in Bravais basis
! SOURCE
!-------------------------------------------------------------------
function ivec_dot_sym(ivec,s)
integer, dimension(3) :: ivec_dot_sym
integer, intent(IN), dimension(:) :: ivec
type(symmetry_type), intent(IN) :: s
!***
if (s%basis == 'Bravais ') then
ivec_dot_sym = nint( ivec_dot_mat(ivec,s%m) )
else
write(6,*) 'Improper use of ivec_dot_sym - not a Bravais symmetry'
stop
endif
end function ivec_dot_sym
!===================================================================
!------------------------------------------------------------------
!****ip group_mod/sym_dot_sym
! FUNCTION
! type(symmetry_type) sym_dot_sym(v,s)
! RETURN
! Dot product of two symmetries
! SOURCE
!-------------------------------------------------------------------
type(symmetry_type) function sym_dot_sym(s1,s2)
type(symmetry_type), intent(IN) :: s1,s2
!***
type(symmetry_type) :: s3
! Variable s3 added to compile on Regatta (6/2004)
if (s1%basis == s2%basis) then
s3%m = mat_dot_mat(s1%m,s2%m)
s3%v = mat_dot_vec(s1%m,s2%v) + s1%v
s3%basis = s1%basis
else
write(6,*) 'Incompatible bases in sym_dot_sym'
stop
endif
sym_dot_sym = s3
end function sym_dot_sym
!===================================================================
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
! Definitions of inverse
!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
!-------------------------------------------------------------------
!****ip group_mod/matrix_inverse
! FUNCTION
! matrix_inverse(m)
! real(long) dimension(3,3) :: matrix_inverse(3,3)
! RETURN
! Matrix inverse for dim=1, dim=2 or dim=3
! SOURCE
!-------------------------------------------------------------------
function matrix_inverse(m)
real(long) :: matrix_inverse(3,3)
real(long), intent(IN):: m(3,3)
!***
real(long) :: cofac(3,3)
real(long):: det
integer :: i, j
matrix_inverse = 0.0_long
if (dim == 3) then
cofac(1,1) = m(2,2)*m(3,3) - m(3,2)*m(2,3)
cofac(1,2) = m(2,1)*m(3,3) - m(3,1)*m(2,3)
cofac(1,3) = m(2,1)*m(3,2) - m(3,1)*m(2,2)
cofac(2,1) = m(1,2)*m(3,3) - m(3,2)*m(1,3)
cofac(2,2) = m(1,1)*m(3,3) - m(3,1)*m(1,3)
cofac(2,3) = m(1,1)*m(3,2) - m(3,1)*m(1,2)
cofac(3,1) = m(1,2)*m(2,3) - m(2,2)*m(1,3)
cofac(3,2) = m(1,1)*m(2,3) - m(2,1)*m(1,3)
cofac(3,3) = m(1,1)*m(2,2) - m(2,1)*m(1,2)
cofac(1,2) = -cofac(1,2)
cofac(2,3) = -cofac(2,3)
cofac(2,1) = -cofac(2,1)
cofac(3,2) = -cofac(3,2)
det = 0.0_long
do i=1, 3
det = det + m(1,i)*cofac(1,i)
enddo
do i=1, 3
do j=1, 3
matrix_inverse(i,j) = cofac(j,i)/det
enddo
enddo
else if (dim == 2) then
cofac(1,1) = m(2,2)
cofac(2,2) = m(1,1)
cofac(1,2) = -m(2,1)
cofac(2,1) = -m(1,2)
det = m(1,1)*m(2,2) - m(1,2)*m(2,1)
do i=1, 2
do j=1, 2
matrix_inverse(i,j) = cofac(j,i)/det
enddo
enddo
else if (dim == 1) then
matrix_inverse(1,1) = 1.0_long/m(1,1)
else
write(6,*) 'Invalid dim in matrix_inverse'
stop
endif
end function matrix_inverse
!===================================================================
!-------------------------------------------------------------------
!****ip group_mod/symmetry_inverse
! FUNCTION
! type(symmetry_type) symmetry_inverse(m)
! RETURN
! Inverse of symmetry operation s
! SOURCE
!-------------------------------------------------------------------
type(symmetry_type) function symmetry_inverse(s)
type(symmetry_type), intent(IN) :: s
!***
real(long) :: inv_m(3,3)
! Inverse stored in variable inv_m to compile on Regatta (6/2004)
symmetry_inverse%basis = s%basis
inv_m = matrix_inverse(s%m)
symmetry_inverse%m = inv_m
symmetry_inverse%v = (-1.0d0)*(symmetry_inverse%m .dot. s%v)
end function symmetry_inverse
!===================================================================