TABLE OF CONTENTS


ABINIT/m_spgbuilder [ Modules ]

[ Top ] [ Modules ]

NAME

 m_spgbuilder

FUNCTION

  Spacegroup builder.

COPYRIGHT

  Copyright (C) 2008-2024 ABINIT group (RC, XG)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_spgbuilder
23 
24  use defs_basis
25  use m_errors
26  use m_abicore
27 
28  use m_symtk,   only : sg_multable, print_symmetries
29  use m_symsg,   only : symsgcube, symsghexa, symsgmono, symsgortho, symsgtetra
30 
31  implicit none
32 
33  private

m_spgbuilder/gensymshub [ Functions ]

[ Top ] [ m_spgbuilder ] [ Functions ]

NAME

 gensymshub

FUNCTION

 Analyse the Shubnikov space group, from the input of the
 the Fedorov space group number, and the Shubnikov space group number:
 1) determine the type (III or IV)
 2) for type (IV), determine the translation generating
   the anti-ferromagnetic operations
 At present, follow strictly the specifications of Bradley
 and Cracknell. However, one should take into account the
 orientation of the symmetry group (spgaxor).

INPUTS

 spgroup = number of space group
 spgroupma = number of magnetic space group

OUTPUT

 genafm(3) = in case of shubnikov type IV, translation, generator of the
  anti-ferromagnetic symmetry operations
 shubnikov = type of the shubnikov group

SOURCE

 580 subroutine gensymshub(genafm,spgroup,spgroupma,shubnikov)
 581 
 582 !Arguments ------------------------------------
 583 !scalars
 584  integer,intent(in) :: spgroup,spgroupma
 585  integer,intent(out) :: shubnikov
 586 !arrays
 587  real(dp),intent(out) :: genafm(3)
 588 
 589 !Local variables ------------------------------
 590 !scalars
 591  integer :: brvlttbw=0,spgrmatch=1
 592  character(len=500) :: message
 593 
 594 ! *************************************************************************
 595 
 596 !List of the input parameters
 597 !DEBUG
 598 !write(std_out,*)
 599 !write(std_out,*)' gensymshub : enter with:'
 600 !write(std_out,*)' brvlttbw  = ',brvlttbw
 601 !write(std_out,*)' spgroup = ',spgroup
 602 !write(std_out,*)' spgroupma = ',spgroupma
 603 !ENDDEBUG
 604 
 605 !Test for consistency the magnetic and non-magnetic space group
 606  select case (spgroup)
 607  case(1)
 608    if (.not.(spgroupma>=3  .and.  3 >=spgroupma) ) spgrmatch=0
 609  case(2)
 610    if (.not.(spgroupma>=6  .and.  7 >=spgroupma) ) spgrmatch=0
 611  case(3)
 612    if (.not.(spgroupma>=3   .and.  6 >=spgroupma) ) spgrmatch=0
 613  case(4)
 614    if (.not.(spgroupma>=9   .and.  12 >=spgroupma) ) spgrmatch=0
 615  case(5)
 616    if (.not.(spgroupma>=15   .and.  17 >=spgroupma) ) spgrmatch=0
 617  case(6)
 618    if (.not.(spgroupma>=20   .and.  23 >=spgroupma) ) spgrmatch=0
 619  case(7)
 620    if (.not.(spgroupma>=26   .and.  31 >=spgroupma) ) spgrmatch=0
 621  case(8)
 622    if (.not.(spgroupma>=34   .and.  36 >=spgroupma) ) spgrmatch=0
 623  case(9)
 624    if (.not.(spgroupma>=39   .and.  41 >=spgroupma) ) spgrmatch=0
 625  case(10)
 626    if (.not.(spgroupma>=44   .and.  49 >=spgroupma) ) spgrmatch=0
 627  case(11)
 628    if (.not.(spgroupma>=52   .and.  57 >=spgroupma) ) spgrmatch=0
 629  case(12)
 630    if (.not.(spgroupma>=60   .and.  64 >=spgroupma) ) spgrmatch=0
 631  case(13)
 632    if (.not.(spgroupma>=67   .and.  74 >=spgroupma) ) spgrmatch=0
 633  case(14)
 634    if (.not.(spgroupma>=77   .and.  84 >=spgroupma) ) spgrmatch=0
 635  case(15)
 636    if (.not.(spgroupma>=87   .and.  91 >=spgroupma) ) spgrmatch=0
 637  case(16)
 638    if (.not.(spgroupma>= 3  .and.   6>=spgroupma) ) spgrmatch=0
 639  case(17)
 640    if (.not.(spgroupma>= 9  .and.  15 >=spgroupma) ) spgrmatch=0
 641  case(18)
 642    if (.not.(spgroupma>=18   .and.  24 >=spgroupma) ) spgrmatch=0
 643  case(19)
 644    if (.not.(spgroupma>=27   .and.  30 >=spgroupma) ) spgrmatch=0
 645  case(20)
 646    if (.not.(spgroupma>=33   .and.  37 >=spgroupma) ) spgrmatch=0
 647  case(21)
 648    if (.not.(spgroupma>=40   .and.  44 >=spgroupma) ) spgrmatch=0
 649  case(22)
 650    if (.not.(spgroupma>=47   .and.  48 >=spgroupma) ) spgrmatch=0
 651  case(23)
 652    if (.not.(spgroupma>=51   .and.  52 >=spgroupma) ) spgrmatch=0
 653  case(24)
 654    if (.not.(spgroupma>=55   .and.  56 >=spgroupma) ) spgrmatch=0
 655  case(25)
 656    if (.not.(spgroupma>=59   .and.  65 >=spgroupma) ) spgrmatch=0
 657  case(26)
 658    if (.not.(spgroupma>=68   .and.  77 >=spgroupma) ) spgrmatch=0
 659  case(27)
 660    if (.not.(spgroupma>=80   .and.  86 >=spgroupma) ) spgrmatch=0
 661  case(28)
 662    if (.not.(spgroupma>=89   .and.  98 >=spgroupma) ) spgrmatch=0
 663  case(29)
 664    if (.not.(spgroupma>=101  .and.  110 >=spgroupma) ) spgrmatch=0
 665  case(30)
 666    if (.not.(spgroupma>=113   .and. 122  >=spgroupma) ) spgrmatch=0
 667  case(31)
 668    if (.not.(spgroupma>=125   .and. 134  >=spgroupma) ) spgrmatch=0
 669  case(32)
 670    if (.not.(spgroupma>=137   .and. 143  >=spgroupma) ) spgrmatch=0
 671  case(33)
 672    if (.not.(spgroupma>=146   .and.  155 >=spgroupma) ) spgrmatch=0
 673  case(34)
 674    if (.not.(spgroupma>=158   .and.  164 >=spgroupma) ) spgrmatch=0
 675  case(35)
 676    if (.not.(spgroupma>=167   .and.   171>=spgroupma) ) spgrmatch=0
 677  case(36)
 678    if (.not.(spgroupma>=174   .and.   179>=spgroupma) ) spgrmatch=0
 679  case(37)
 680    if (.not.(spgroupma>=182   .and.  186 >=spgroupma) ) spgrmatch=0
 681  case(38)
 682    if (.not.(spgroupma>=189   .and.  194 >=spgroupma) ) spgrmatch=0
 683  case(39)
 684    if (.not.(spgroupma>=197   .and.  202 >=spgroupma) ) spgrmatch=0
 685  case(40)
 686    if (.not.(spgroupma>=205   .and.  210 >=spgroupma) ) spgrmatch=0
 687  case(41)
 688    if (.not.(spgroupma>=213   .and.  218 >=spgroupma) ) spgrmatch=0
 689  case(42)
 690    if (.not.(spgroupma>=221   .and.   223>=spgroupma) ) spgrmatch=0
 691  case(43)
 692    if (.not.(spgroupma>=226   .and.   228>=spgroupma) ) spgrmatch=0
 693  case(44)
 694    if (.not.(spgroupma>=231   .and.  234 >=spgroupma) ) spgrmatch=0
 695  case(45)
 696    if (.not.(spgroupma>=237   .and.  240 >=spgroupma) ) spgrmatch=0
 697  case(46)
 698    if (.not.(spgroupma>=243   .and.  248 >=spgroupma) ) spgrmatch=0
 699  case(47)
 700    if (.not.(spgroupma>=251   .and.  256 >=spgroupma) ) spgrmatch=0
 701  case(48)
 702    if (.not.(spgroupma>=259   .and.  264 >=spgroupma) ) spgrmatch=0
 703  case(49)
 704    if (.not.(spgroupma>=267   .and.  276 >=spgroupma) ) spgrmatch=0
 705  case(50)
 706    if (.not.(spgroupma>=279   .and.  288 >=spgroupma) ) spgrmatch=0
 707  case(51)
 708    if (.not.(spgroupma>=291   .and.  304 >=spgroupma) ) spgrmatch=0
 709  case(52)
 710    if (.not.(spgroupma>=307   .and.  320 >=spgroupma) ) spgrmatch=0
 711  case(53)
 712    if (.not.(spgroupma>=323   .and.  336 >=spgroupma) ) spgrmatch=0
 713  case(54)
 714    if (.not.(spgroupma>=339   .and.  352 >=spgroupma) ) spgrmatch=0
 715  case(55)
 716    if (.not.(spgroupma>=355   .and.  364 >=spgroupma) ) spgrmatch=0
 717  case(56)
 718    if (.not.(spgroupma>=367   .and.  376 >=spgroupma) ) spgrmatch=0
 719  case(57)
 720    if (.not.(spgroupma>=379   .and.  392 >=spgroupma) ) spgrmatch=0
 721  case(58)
 722    if (.not.(spgroupma>=395   .and.  404 >=spgroupma) ) spgrmatch=0
 723  case(59)
 724    if (.not.(spgroupma>=407   .and.  416 >=spgroupma) ) spgrmatch=0
 725  case(60)
 726    if (.not.(spgroupma>=419   .and.  432 >=spgroupma) ) spgrmatch=0
 727  case(61)
 728    if (.not.(spgroupma>=435   .and.  440 >=spgroupma) ) spgrmatch=0
 729  case(62)
 730    if (.not.(spgroupma>=443   .and.  456 >=spgroupma) ) spgrmatch=0
 731  case(63)
 732    if (.not.(spgroupma>=459   .and.  468 >=spgroupma) ) spgrmatch=0
 733  case(64)
 734    if (.not.(spgroupma>=471   .and.  480 >=spgroupma) ) spgrmatch=0
 735  case(65)
 736    if (.not.(spgroupma>=483   .and.  490 >=spgroupma) ) spgrmatch=0
 737  case(66)
 738    if (.not.(spgroupma>=493   .and.  500 >=spgroupma) ) spgrmatch=0
 739  case(67)
 740    if (.not.(spgroupma>=503   .and.  510 >=spgroupma) ) spgrmatch=0
 741  case(68)
 742    if (.not.(spgroupma>=513   .and.  520 >=spgroupma) ) spgrmatch=0
 743  case(69)
 744    if (.not.(spgroupma>=523   .and.  526 >=spgroupma) ) spgrmatch=0
 745  case(70)
 746    if (.not.(spgroupma>=529   .and.  532 >=spgroupma) ) spgrmatch=0
 747  case(71)
 748    if (.not.(spgroupma>=535   .and.  538 >=spgroupma) ) spgrmatch=0
 749  case(72)
 750    if (.not.(spgroupma>=541   .and.  547 >=spgroupma) ) spgrmatch=0
 751  case(73)
 752    if (.not.(spgroupma>=550   .and.  553 >=spgroupma) ) spgrmatch=0
 753  case(74)
 754    if (.not.(spgroupma>=556   .and.  562 >=spgroupma) ) spgrmatch=0
 755  case(75)
 756    if (.not.(spgroupma>= 3  .and.   6>=spgroupma) ) spgrmatch=0
 757  case(76)
 758    if (.not.(spgroupma>= 9  .and.  12 >=spgroupma) ) spgrmatch=0
 759  case(77)
 760    if (.not.(spgroupma>= 15  .and.   18>=spgroupma) ) spgrmatch=0
 761  case(78)
 762    if (.not.(spgroupma>= 21  .and.   24>=spgroupma) ) spgrmatch=0
 763  case(79)
 764    if (.not.(spgroupma>= 27  .and.  28 >=spgroupma) ) spgrmatch=0
 765  case(80)
 766    if (.not.(spgroupma>= 31  .and.  32 >=spgroupma) ) spgrmatch=0
 767  case(81)
 768    if (.not.(spgroupma>= 35  .and.  38 >=spgroupma) ) spgrmatch=0
 769  case(82)
 770    if (.not.(spgroupma>= 41  .and.  42 >=spgroupma) ) spgrmatch=0
 771  case(83)
 772    if (.not.(spgroupma>=45   .and.   50>=spgroupma) ) spgrmatch=0
 773  case(84)
 774    if (.not.(spgroupma>=53   .and.   58>=spgroupma) ) spgrmatch=0
 775  case(85)
 776    if (.not.(spgroupma>=61   .and.   66>=spgroupma) ) spgrmatch=0
 777  case(86)
 778    if (.not.(spgroupma>=69   .and.   74>=spgroupma) ) spgrmatch=0
 779  case(87)
 780    if (.not.(spgroupma>=77   .and.   80>=spgroupma) ) spgrmatch=0
 781  case(88)
 782    if (.not.(spgroupma>=83   .and.   86>=spgroupma) ) spgrmatch=0
 783  case(89)
 784    if (.not.(spgroupma>=89   .and.   94>=spgroupma) ) spgrmatch=0
 785  case(90)
 786    if (.not.(spgroupma>=97   .and.   102>=spgroupma) ) spgrmatch=0
 787  case(91)
 788    if (.not.(spgroupma>=105   .and.   110>=spgroupma) ) spgrmatch=0
 789  case(92)
 790    if (.not.(spgroupma>=113   .and.   118>=spgroupma) ) spgrmatch=0
 791  case(93)
 792    if (.not.(spgroupma>=121   .and.   126>=spgroupma) ) spgrmatch=0
 793  case(94)
 794    if (.not.(spgroupma>=129   .and.   134>=spgroupma) ) spgrmatch=0
 795  case(95)
 796    if (.not.(spgroupma>=137   .and.   142>=spgroupma) ) spgrmatch=0
 797  case(96)
 798    if (.not.(spgroupma>=145   .and.   150>=spgroupma) ) spgrmatch=0
 799  case(97)
 800    if (.not.(spgroupma>=153   .and.   156>=spgroupma) ) spgrmatch=0
 801  case(98)
 802    if (.not.(spgroupma>=159   .and.   162>=spgroupma) ) spgrmatch=0
 803  case(99)
 804    if (.not.(spgroupma>=165   .and.   170>=spgroupma) ) spgrmatch=0
 805  case(100)
 806    if (.not.(spgroupma>=173   .and.   178>=spgroupma) ) spgrmatch=0
 807  case(101)
 808    if (.not.(spgroupma>=181   .and.   186>=spgroupma) ) spgrmatch=0
 809  case(102)
 810    if (.not.(spgroupma>=189   .and.   194>=spgroupma) ) spgrmatch=0
 811  case(103)
 812    if (.not.(spgroupma>=197   .and.   202>=spgroupma) ) spgrmatch=0
 813  case(104)
 814    if (.not.(spgroupma>=205   .and.   210>=spgroupma) ) spgrmatch=0
 815  case(105)
 816    if (.not.(spgroupma>=213   .and.   218>=spgroupma) ) spgrmatch=0
 817  case(106)
 818    if (.not.(spgroupma>=221   .and.   226>=spgroupma) ) spgrmatch=0
 819  case(107)
 820    if (.not.(spgroupma>=229  .and.   232>=spgroupma) ) spgrmatch=0
 821  case(108)
 822    if (.not.(spgroupma>=235   .and.   238>=spgroupma) ) spgrmatch=0
 823  case(109)
 824    if (.not.(spgroupma>=241   .and.   244>=spgroupma) ) spgrmatch=0
 825  case(110)
 826    if (.not.(spgroupma>=247   .and.   250>=spgroupma) ) spgrmatch=0
 827  case(111)
 828    if (.not.(spgroupma>=253   .and.   258>=spgroupma) ) spgrmatch=0
 829  case(112)
 830    if (.not.(spgroupma>=261   .and.   266>=spgroupma) ) spgrmatch=0
 831  case(113)
 832    if (.not.(spgroupma>=269   .and.   274>=spgroupma) ) spgrmatch=0
 833  case(114)
 834    if (.not.(spgroupma>=277   .and.   282>=spgroupma) ) spgrmatch=0
 835  case(115)
 836    if (.not.(spgroupma>=285   .and.   290>=spgroupma) ) spgrmatch=0
 837  case(116)
 838    if (.not.(spgroupma>=293   .and.   298>=spgroupma) ) spgrmatch=0
 839  case(117)
 840    if (.not.(spgroupma>=301   .and.   306>=spgroupma) ) spgrmatch=0
 841  case(118)
 842    if (.not.(spgroupma>=309   .and.   314>=spgroupma) ) spgrmatch=0
 843  case(119)
 844    if (.not.(spgroupma>=317   .and.   320>=spgroupma) ) spgrmatch=0
 845  case(120)
 846    if (.not.(spgroupma>=323   .and.   326>=spgroupma) ) spgrmatch=0
 847  case(121)
 848    if (.not.(spgroupma>=329   .and.   332>=spgroupma) ) spgrmatch=0
 849  case(122)
 850    if (.not.(spgroupma>=335   .and.   338>=spgroupma) ) spgrmatch=0
 851  case(123)
 852    if (.not.(spgroupma>=341   .and.   350>=spgroupma) ) spgrmatch=0
 853  case(124)
 854    if (.not.(spgroupma>=353   .and.   362>=spgroupma) ) spgrmatch=0
 855  case(125)
 856    if (.not.(spgroupma>=365   .and.   374>=spgroupma) ) spgrmatch=0
 857  case(126)
 858    if (.not.(spgroupma>=377   .and.   386>=spgroupma) ) spgrmatch=0
 859  case(127)
 860    if (.not.(spgroupma>=389   .and.   398>=spgroupma) ) spgrmatch=0
 861  case(128)
 862    if (.not.(spgroupma>=401   .and.   410>=spgroupma) ) spgrmatch=0
 863  case(129)
 864    if (.not.(spgroupma>=413   .and.   422>=spgroupma) ) spgrmatch=0
 865  case(130)
 866    if (.not.(spgroupma>=425   .and.   434>=spgroupma) ) spgrmatch=0
 867  case(131)
 868    if (.not.(spgroupma>=437   .and.   446>=spgroupma) ) spgrmatch=0
 869  case(132)
 870    if (.not.(spgroupma>=449   .and.   458>=spgroupma) ) spgrmatch=0
 871  case(133)
 872    if (.not.(spgroupma>=461   .and.   470>=spgroupma) ) spgrmatch=0
 873  case(134)
 874    if (.not.(spgroupma>=473   .and.   482>=spgroupma) ) spgrmatch=0
 875  case(135)
 876    if (.not.(spgroupma>=485   .and.   494>=spgroupma) ) spgrmatch=0
 877  case(136)
 878    if (.not.(spgroupma>=497  .and.   506>=spgroupma) ) spgrmatch=0
 879  case(137)
 880    if (.not.(spgroupma>=509   .and.   518>=spgroupma) ) spgrmatch=0
 881  case(138)
 882    if (.not.(spgroupma>=521   .and.   530>=spgroupma) ) spgrmatch=0
 883  case(139)
 884    if (.not.(spgroupma>=533   .and.   540>=spgroupma) ) spgrmatch=0
 885  case(140)
 886    if (.not.(spgroupma>=543   .and.   550>=spgroupma) ) spgrmatch=0
 887  case(141)
 888    if (.not.(spgroupma>=553   .and.   560>=spgroupma) ) spgrmatch=0
 889  case(142)
 890    if (.not.(spgroupma>=563   .and.   570>=spgroupma) ) spgrmatch=0
 891  case(143)
 892    if (.not.(spgroupma>=3   .and.   3>=spgroupma) ) spgrmatch=0
 893  case(144)
 894    if (.not.(spgroupma>=6   .and.   6>=spgroupma) ) spgrmatch=0
 895  case(145)
 896    if (.not.(spgroupma>=9   .and.   9>=spgroupma) ) spgrmatch=0
 897  case(146)
 898    if (.not.(spgroupma>=12   .and.   12>=spgroupma) ) spgrmatch=0
 899  case(147)
 900    if (.not.(spgroupma>=15   .and.   16>=spgroupma) ) spgrmatch=0
 901  case(148)
 902    if (.not.(spgroupma>=19   .and.   20>=spgroupma) ) spgrmatch=0
 903  case(149)
 904    if (.not.(spgroupma>=23   .and.   24>=spgroupma) ) spgrmatch=0
 905  case(150)
 906    if (.not.(spgroupma>=27   .and.   28>=spgroupma) ) spgrmatch=0
 907  case(151)
 908    if (.not.(spgroupma>=31   .and.   32>=spgroupma) ) spgrmatch=0
 909  case(152)
 910    if (.not.(spgroupma>=35   .and.   36>=spgroupma) ) spgrmatch=0
 911  case(153)
 912    if (.not.(spgroupma>=39   .and.   40>=spgroupma) ) spgrmatch=0
 913  case(154)
 914    if (.not.(spgroupma>=43   .and.   44>=spgroupma) ) spgrmatch=0
 915  case(155)
 916    if (.not.(spgroupma>=47   .and.   48>=spgroupma) ) spgrmatch=0
 917  case(156)
 918    if (.not.(spgroupma>=51   .and.   52>=spgroupma) ) spgrmatch=0
 919  case(157)
 920    if (.not.(spgroupma>=55   .and.   56>=spgroupma) ) spgrmatch=0
 921  case(158)
 922    if (.not.(spgroupma>=59   .and.   60>=spgroupma) ) spgrmatch=0
 923  case(159)
 924    if (.not.(spgroupma>=63   .and.   64>=spgroupma) ) spgrmatch=0
 925  case(160)
 926    if (.not.(spgroupma>=67   .and.   68>=spgroupma) ) spgrmatch=0
 927  case(161)
 928    if (.not.(spgroupma>=71   .and.   72>=spgroupma) ) spgrmatch=0
 929  case(162)
 930    if (.not.(spgroupma>=75   .and.   78>=spgroupma) ) spgrmatch=0
 931  case(163)
 932    if (.not.(spgroupma>=81   .and.   84>=spgroupma) ) spgrmatch=0
 933  case(164)
 934    if (.not.(spgroupma>=87   .and.   90>=spgroupma) ) spgrmatch=0
 935  case(165)
 936    if (.not.(spgroupma>=93   .and.   96>=spgroupma) ) spgrmatch=0
 937  case(166)
 938    if (.not.(spgroupma>=99   .and.   102>=spgroupma) ) spgrmatch=0
 939  case(167)
 940    if (.not.(spgroupma>=105   .and.   108>=spgroupma) ) spgrmatch=0
 941  case(168)
 942    if (.not.(spgroupma>=111   .and.   112>=spgroupma) ) spgrmatch=0
 943  case(169)
 944    if (.not.(spgroupma>=115   .and.   116>=spgroupma) ) spgrmatch=0
 945  case(170)
 946    if (.not.(spgroupma>=119   .and.   120>=spgroupma) ) spgrmatch=0
 947  case(171)
 948    if (.not.(spgroupma>=123   .and.   124>=spgroupma) ) spgrmatch=0
 949  case(172)
 950    if (.not.(spgroupma>=127   .and.   128>=spgroupma) ) spgrmatch=0
 951  case(173)
 952    if (.not.(spgroupma>=131   .and.   132>=spgroupma) ) spgrmatch=0
 953  case(174)
 954    if (.not.(spgroupma>=135   .and.   136>=spgroupma) ) spgrmatch=0
 955  case(175)
 956    if (.not.(spgroupma>=139   .and.   142>=spgroupma) ) spgrmatch=0
 957  case(176)
 958    if (.not.(spgroupma>=145   .and.   148>=spgroupma) ) spgrmatch=0
 959  case(177)
 960    if (.not.(spgroupma>=151   .and.   154>=spgroupma) ) spgrmatch=0
 961  case(178)
 962    if (.not.(spgroupma>=157   .and.   160>=spgroupma) ) spgrmatch=0
 963  case(179)
 964    if (.not.(spgroupma>=163   .and.   166>=spgroupma) ) spgrmatch=0
 965  case(180)
 966    if (.not.(spgroupma>=169   .and.   172>=spgroupma) ) spgrmatch=0
 967  case(181)
 968    if (.not.(spgroupma>=175   .and.   178>=spgroupma) ) spgrmatch=0
 969  case(182)
 970    if (.not.(spgroupma>=181   .and.   184>=spgroupma) ) spgrmatch=0
 971  case(183)
 972    if (.not.(spgroupma>=187   .and.   190>=spgroupma) ) spgrmatch=0
 973  case(184)
 974    if (.not.(spgroupma>=193   .and.   196>=spgroupma) ) spgrmatch=0
 975  case(185)
 976    if (.not.(spgroupma>=199   .and.   202>=spgroupma) ) spgrmatch=0
 977  case(186)
 978    if (.not.(spgroupma>=205   .and.   208>=spgroupma) ) spgrmatch=0
 979  case(187)
 980    if (.not.(spgroupma>=211   .and.   214>=spgroupma) ) spgrmatch=0
 981  case(188)
 982    if (.not.(spgroupma>=217   .and.   220>=spgroupma) ) spgrmatch=0
 983  case(189)
 984    if (.not.(spgroupma>=223   .and.   226>=spgroupma) ) spgrmatch=0
 985  case(190)
 986    if (.not.(spgroupma>=229   .and.   232>=spgroupma) ) spgrmatch=0
 987  case(191)
 988    if (.not.(spgroupma>=235   .and.   242>=spgroupma) ) spgrmatch=0
 989  case(192)
 990    if (.not.(spgroupma>=245   .and.   252>=spgroupma) ) spgrmatch=0
 991  case(193)
 992    if (.not.(spgroupma>=255   .and.   262>=spgroupma) ) spgrmatch=0
 993  case(194)
 994    if (.not.(spgroupma>=265   .and.   272>=spgroupma) ) spgrmatch=0
 995  case(195)
 996    if (.not.(spgroupma>=3   .and.   3>=spgroupma) ) spgrmatch=0
 997  case(196)
 998    if (.not.(spgroupma>=6   .and.   6>=spgroupma) ) spgrmatch=0
 999  case(197)
1000    spgrmatch=0
1001  case(198)
1002    if (.not.(spgroupma>=11  .and.   11>=spgroupma) ) spgrmatch=0
1003  case(199)
1004    spgrmatch=0
1005  case(200)
1006    if (.not.(spgroupma>=16   .and.   17>=spgroupma) ) spgrmatch=0
1007  case(201)
1008    if (.not.(spgroupma>=20   .and.   21>=spgroupma) ) spgrmatch=0
1009  case(202)
1010    if (.not.(spgroupma>=24   .and.   25>=spgroupma) ) spgrmatch=0
1011  case(203)
1012    if (.not.(spgroupma>=28   .and.   29>=spgroupma) ) spgrmatch=0
1013  case(204)
1014    if (.not.(spgroupma>=32   .and.   32>=spgroupma) ) spgrmatch=0
1015  case(205)
1016    if (.not.(spgroupma>=35   .and.   36>=spgroupma) ) spgrmatch=0
1017  case(206)
1018    if (.not.(spgroupma>=39   .and.   39>=spgroupma) ) spgrmatch=0
1019  case(207)
1020    if (.not.(spgroupma>=42   .and.   43>=spgroupma) ) spgrmatch=0
1021  case(208)
1022    if (.not.(spgroupma>=46   .and.   47>=spgroupma) ) spgrmatch=0
1023  case(209)
1024    if (.not.(spgroupma>=50   .and.   51>=spgroupma) ) spgrmatch=0
1025  case(210)
1026    if (.not.(spgroupma>=54   .and.   55>=spgroupma) ) spgrmatch=0
1027  case(211)
1028    if (.not.(spgroupma>=58   .and.   58>=spgroupma) ) spgrmatch=0
1029  case(212)
1030    if (.not.(spgroupma>=61  .and.   62>=spgroupma) ) spgrmatch=0
1031  case(213)
1032    if (.not.(spgroupma>=65   .and.   66>=spgroupma) ) spgrmatch=0
1033  case(214)
1034    if (.not.(spgroupma>=69   .and.   69>=spgroupma) ) spgrmatch=0
1035  case(215)
1036    if (.not.(spgroupma>=72   .and.   73>=spgroupma) ) spgrmatch=0
1037  case(216)
1038    if (.not.(spgroupma>=76   .and.   77>=spgroupma) ) spgrmatch=0
1039  case(217)
1040    if (.not.(spgroupma>=80   .and.   80>=spgroupma) ) spgrmatch=0
1041  case(218)
1042    if (.not.(spgroupma>=83   .and.   84>=spgroupma) ) spgrmatch=0
1043  case(219)
1044    if (.not.(spgroupma>=87   .and.   88>=spgroupma) ) spgrmatch=0
1045  case(220)
1046    if (.not.(spgroupma>=91   .and.   91>=spgroupma) ) spgrmatch=0
1047  case(221)
1048    if (.not.(spgroupma>=94   .and.   97>=spgroupma) ) spgrmatch=0
1049  case(222)
1050    if (.not.(spgroupma>=100  .and.   103>=spgroupma) ) spgrmatch=0
1051  case(223)
1052    if (.not.(spgroupma>=106   .and.   109>=spgroupma) ) spgrmatch=0
1053  case(224)
1054    if (.not.(spgroupma>=112   .and.   115>=spgroupma) ) spgrmatch=0
1055  case(225)
1056    if (.not.(spgroupma>=118   .and.   121>=spgroupma) ) spgrmatch=0
1057  case(226)
1058    if (.not.(spgroupma>=124   .and.   127>=spgroupma) ) spgrmatch=0
1059  case(227)
1060    if (.not.(spgroupma>=130   .and.   133>=spgroupma) ) spgrmatch=0
1061  case(228)
1062    if (.not.(spgroupma>=136   .and.   139>=spgroupma) ) spgrmatch=0
1063  case(229)
1064    if (.not.(spgroupma>=142   .and.   144>=spgroupma) ) spgrmatch=0
1065  case(230)
1066    if (.not.(spgroupma>=147   .and.   149>=spgroupma) ) spgrmatch=0
1067  case default
1068    write(message, '(3a,i8,4a)' )&
1069 &   'The non-magnetic spacegroup is not specified ',ch10,&
1070 &   'while the magnetic space group is specified, spgroupma= ',spgroupma,ch10,&
1071 &   'This is not allowed.  ',ch10,&
1072 &   'Action: specify spgroup in the input file.'
1073    ABI_ERROR(message)
1074  end select
1075 
1076  if (spgrmatch==0) then
1077    write(message, '(a,i8,a,a,i8,4a)' )&
1078 &   'mismatch between the non-magnetic spacegroup ',spgroup,ch10,&
1079 &   'and the magnetic space group ',spgroupma,ch10,&
1080 &   'This is not allowed.  ',ch10,&
1081 &   'Action: modify spgroup or spgroupma in the input file.'
1082    ABI_ERROR(message)
1083  end if
1084 
1085 !DEBUG
1086 !write(std_out,*) ' gensymshub, after check the spgroup ... '
1087 !ENDDEBUG
1088 
1089 !Assign the magnetic Bravais lattice type from the magnetic space group number
1090 !As the magnetic space group number begins with 1 for EACH crystal system
1091 !we must first make our choice as a function of spgroup
1092 
1093 !Convention :
1094 !brvlttbw = input variable giving Bravais black-and-white translation
1095 !(from 1 to 8 : Shubnikov type IV space group)
1096 !1,2,3 = 1/2 translation along a, b, or c respectively
1097 !4,5,6 = translation corresponding to A,B,C centering, respectively
1098 !7 = I centering
1099 !8 = (1/2 0 0) centering corresponding to normal F lattice
1100 !9 -> Shubnikov type III space group
1101 !Note that the use of the table 7.3 (p585) of Bradney and Cracknell
1102 !for the definition of the translation vectors
1103 !is extremely confusing, due to the strange choice of basis
1104 !vectors of table 3.1.
1105 !See table 7.4 (p588) for the spgroupma interpretation
1106 
1107  select case(spgroup)
1108  case(1,2)     !Triclinic
1109    select case(spgroupma)
1110    case(6)
1111      brvlttbw=9      !ShubIII
1112    case(3,7)
1113      brvlttbw=1      !Ps (note that it is not body centered, according to Table
1114 !        7.3 of Bradley and Cracknell)
1115    end select
1116  case(3:15)     !Monoclinic
1117    select case(spgroupma)
1118    case(3,9,15,20,26,34,39,44:46,52:54,60:62,67:69,77:79,87:89)
1119      brvlttbw=9      !ShubIII
1120    case(4,10,17,21,27,36,41,47,55,64,70,80,91)
1121      brvlttbw=1     !a
1122    case(5,11,22,29,48,56,71,81)
1123      brvlttbw=2     !b
1124    case(16,28,35,40,63,72,82,90)
1125      brvlttbw=3     !c
1126    case(31,73,83)
1127      brvlttbw=4     !A
1128    case(6,12,23,30,49,57,74,84)
1129      brvlttbw=6     !C
1130    end select
1131  case(16:74)     !Orthorhombic
1132    select case(spgroupma)
1133    case(3,9,10,18,19,27,33,34,40,41,47,51,55,59,60,68:70,&
1134 &     80,81,89:91,101:103,113:115,125:127,137,138,146:148,158,159,&
1135 &     167,168,174:176,182,183,189:191,197:199,205:207,213:215,221,&
1136 &     222,226,227,231,232,237,238,243:245,251:253,259:261,267:271,&
1137 &     279:283,291:297,307:313,323:329,339:345,355:359,367:371,&
1138 &     379:385,395:399,407:411,419:425,435:437,443:449,459:465,&
1139 &     471:477,483:487,493:497,503:507,513:517,523:525,529:531,&
1140 &     535:537,541:545,550:552,556:560)
1141      brvlttbw=9      !ShubIII
1142    case(4,11,20,28,36,43,62,71,83,92,104,116,128,140,&
1143 &     149,160,170,178,185,192,200,208,216,234,240,247,254,262,272,&
1144 &     284,298,314,330,346,360,372,386,400,412,426,438,450,467,&
1145 &     479,489,499,509,519,547,562)
1146      brvlttbw=1     !a
1147    case(72,93,105,117,129,150,248,299,315,331,347,387,427,451)
1148      brvlttbw=2     !b
1149    case(12,21,35,42,52,56,61,73,82,94,106,118,130,139,&
1150 &     151,161,169,177,184,193,201,209,217,233,239,246,273,285,300,&
1151 &     316,332,348,361,373,388,401,413,428,452,466,478,488,498,&
1152 &     508,518,538,546,553,561)
1153      brvlttbw=3     !c
1154    case(13,22,37,44,64,74,85,95,107,119,131,142,152,162,&
1155 &     171,179,186,274,286,301,317,333,349,362,374,389,402,414,429,&
1156 &     453,468,480,490,500,510,520)
1157      brvlttbw=4     !A
1158    case(75,96,108,120,132,153,302,318,334,350,390,430,454)
1159      brvlttbw=5     !B
1160    case(5,14,23,29,63,76,84,97,109,121,133,141,154,163,&
1161 &     194,202,210,218,255,263,275,287,303,319,335,351,363,375,&
1162 &     391,403,415,431,439,455)
1163      brvlttbw=6     !C
1164    case(6,15,24,30,65,77,86,98,110,122,134,143,155,164,256,&
1165 &     264,276,288,304,320,336,352,364,376,392,404,416,432,440,456)
1166      brvlttbw=7     !I
1167    case(48,223,228,526,532)
1168      brvlttbw=8     !s
1169    end select
1170  case(75:142)       !Tetragonal
1171    select case(spgroupma)
1172    case(3,9,15,21,27,31,35,41,45:47,53:55,61:63,69:71,&
1173 &     77:79,83:85,89:91,97:99,105:107,113:115,121:123,129:131,&
1174 &     137:139,145:147,153:155,159:161,165:167,173:175,181:183,&
1175 &     189:191,197:199,205:207,213:215,221:223,229:231,235:237,&
1176 &     241:243,247:249,253:255,261:263,269:271,277:279,285:287,&
1177 &     293:295,301:303,309:311,317:319,323:325,329:331,335:337,&
1178 &     341:347,353:359,365:371,377:383,389:395,401:407,413:419,&
1179 &     425:431,437:443,449:455,461:467,473:479,485:491,497:503,&
1180 &     509:515,521:527,533:539,543:549,553:559,563:569)
1181      brvlttbw=9      !ShubIII
1182    case(4,10,16,22,28,32,36,42,48,56,64,72,80,86,92,100,&
1183 &     108,116,124,132,140,148,156,162,168,176,184,192,200,208,&
1184 &     216,224,232,238,244,250,256,264,272,280,288,296,304,312,&
1185 &     320,326,332,338,348,360,372,384,396,408,420,432,444,456,&
1186 &     468,480,492,504,516,528,540,550,560,570)
1187      brvlttbw=3     !c
1188    case(5,11,17,23,37,49,57,65,73,93,101,109,117,125,133,141,&
1189 &     149,169,177,185,193,201,209,217,225,257,265,273,281,289,297,&
1190 &     305,313,349,361,373,385,397,409,421,433,445,457,469,481,493,&
1191 &     505,517,529)
1192      brvlttbw=6     !C
1193    case(6,12,18,24,38,50,58,66,74,94,102,110,118,126,134,142,&
1194 &     150,170,178,186,194,202,210,218,226,258,266,274,282,290,298,&
1195 &     306,314,350,362,374,386,398,410,422,434,446,458,470,482,494,&
1196 &     506,518,530)
1197      brvlttbw=7     !I
1198    end select
1199  case(143:194)      !Hexagonal
1200    select case(spgroupma)
1201    case(15,19,23,27,31,35,39,43,47,51,55,59,63,67,71,&
1202 &     75:77,81:83,87:89,93:95,99:101,105:107,111,115,119,123,127,&
1203 &     131,135,139:141,145:147,151:153,157:159,163:165,169:171,&
1204 &     175:177,181:183,187:189,193:195,199:201,205:207,211:213,&
1205 &     217:219,223:225,229:231,235:241,245:251,255:261,265:271)
1206      brvlttbw=9      !ShubIII
1207    case(3,6,9,16,24,28,32,36,40,44,52,56,60,64,78,84,&
1208 &     90,96,112,116,120,124,128,132,136,142,148,154,160,166,172,&
1209 &     178,184,190,196,202,208,214,220,226,232,242,252,262,272)
1210 !    brvlttbw=6     !C    XG230719 This is likely erroneous   
1211      brvlttbw=3     !c
1212    case(12,20,48,68,72,102,108)
1213      brvlttbw=7     !I
1214    end select
1215  case(195:230)       !Cubic
1216    select case(spgroupma)
1217    case(16,20,24,28,32,35,39,42,46,50,54,58,61,65,69,&
1218 &     72,76,80,83,87,91,94:96,100:102,106:108,112:114,118:120,&
1219 &     124:126,130:132,136:138,142:144,147:149)
1220      brvlttbw=9      !ShubIII
1221    case(3,11,17,21,36,43,47,62,66,73,84,97,103,109,115)
1222      brvlttbw=7     !I
1223    case(6,25,29,51,55,77,88,121,127,133,139)
1224      brvlttbw=8     !s
1225    end select
1226  end select
1227 
1228  if(brvlttbw==9)shubnikov=3 !Shubnikov type III
1229  if(brvlttbw>=1 .and. brvlttbw<=8) shubnikov=4 !Shubnikov type IV
1230 
1231  genafm(:)=zero
1232  if(shubnikov==4)then
1233    if(brvlttbw==1)genafm(1)=half
1234    if(brvlttbw==2)genafm(2)=half
1235    if(brvlttbw==3)genafm(3)=half
1236    if(brvlttbw>=4 .and. brvlttbw<=8)then
1237      genafm(:)=half
1238      if(brvlttbw==4)genafm(1)=zero
1239      if(brvlttbw==5)genafm(2)=zero
1240      if(brvlttbw==6)genafm(3)=zero
1241    end if
1242  end if
1243 
1244 !DEBUG
1245 !write(std_out,*) 'gensymshub : end '
1246 !write(std_out,*) 'gensymshub : brvlttbw=',brvlttbw
1247 !write(std_out,*) 'gensymshub : genafm=',genafm
1248 !write(std_out,*) 'gensymshub, shubnikov =',shubnikov
1249 !ENDDEBUG
1250 
1251 end subroutine gensymshub

m_spgbuilder/gensymshub4 [ Functions ]

[ Top ] [ m_spgbuilder ] [ Functions ]

NAME

 gensymshub4

FUNCTION

 Assigns the Bravais magnetic translations for shubnikov type IV
 symmetry groups starting from the Fedorov space group
 and the translation, generator of the anti-ferromagnetic
 operations (as input). It will double nsym and tnons. In the end both
 symrel and tnons are completely determined.

INPUTS

 genafm(3) = translation, generator of the anti-ferromagnetic symmatry operations
 msym = maximum number of symmetry operations

OUTPUT

 symafm(msym)= (anti)ferromagnetic part of symmetry operations

SIDE EFFECTS

 nsym=number of symmetry operations, without magnetic operations at input,
  and with magnetic operations at output
 symrel(3,3,msym)=symmetry operations in real space in terms
  of primitive translations, without magnetic operations at input,
  and with magnetic operations at output
 tnons(3,msym)=nonsymmorphic translations for symmetry operations
  without magnetic operations at input,
  and with magnetic operations at output

SOURCE

1284 subroutine gensymshub4(genafm,msym,nsym,symafm,symrel,tnons)
1285 
1286 !Arguments ------------------------------------
1287 !scalars
1288  integer,intent(in) :: msym
1289  integer,intent(inout) :: nsym
1290 !arrays
1291  integer,intent(inout) :: symafm(msym),symrel(3,3,msym)
1292  real(dp),intent(in) :: genafm(3)
1293  real(dp),intent(inout) :: tnons(3,msym)
1294 
1295 !Local variables ------------------------------
1296 !scalars
1297  integer :: ii
1298  character(len=500) :: message
1299 
1300 ! *************************************************************************
1301 
1302  if(msym<2*nsym)then
1303    write(message, '(3a)' )&
1304 &   'The number of symmetries in the Shubnikov type IV space group',ch10,&
1305 &   'is larger than the maximal allowed number of symmetries.'
1306    ABI_ERROR(message)
1307  end if
1308 
1309  do ii=1,nsym
1310    tnons(:,nsym+ii)=tnons(:,ii)+genafm(:)
1311    symrel(:,:,nsym+ii)=symrel(:,:,ii)
1312    symafm(ii)=1
1313    symafm(nsym+ii)=-1
1314  end do
1315  nsym=nsym*2
1316 
1317 end subroutine gensymshub4

m_spgbuilder/gensymspgr [ Functions ]

[ Top ] [ m_spgbuilder ] [ Functions ]

NAME

 gensymspgr

FUNCTION

 Give all the symmetry operations starting from the space group symbol.
 Suppose we are working in a conventional cell
 If brvltt 0 or positive, the pure translations of the
 Bravais lattice are included as generator of the group.
 In brvltt=-1, no pure translation is present, and the
 cell should be changed from conventional to primitive,
 outside of this routine.
 Treat also Shubnikov type III space groups.

INPUTS

 brvltt = input variable giving Bravais lattice
 msym = default number of symmetry operations
 shubnikov= magnetic type of the space group to be generated
 spgaxor = orientation of the cell axis (might be needed)
 spgorig = second choice of origin for certain groups
 (might be needed if nsym==0)
 spgroup = number of space group
 spgroupma= number of the magnetic space group

OUTPUT

 nsym = number of symmetry operations
 symafm(nsym)=(anti)ferromagnetic part of symmetry operations
 symrel(3,3,nsym)=symmetry operations in real space in terms
                  of primitive translations
 tnons(3,nsym)=nonsymmorphic translations for symmetry operations

SIDE EFFECTS

 brvltt = input variable giving Bravais lattice

SOURCE

 80 subroutine gensymspgr(brvltt,msym,nsym,shubnikov,spgaxor,spgorig,spgroup,spgroupma,symafm,symrel,tnons)
 81 
 82 !Arguments ------------------------------------
 83 !scalars
 84  integer,intent(in) :: msym,shubnikov,spgaxor,spgorig,spgroup,spgroupma
 85  integer,intent(inout) :: brvltt
 86  integer,intent(out) :: nsym
 87 !arrays
 88  integer,intent(inout) :: symafm(msym),symrel(3,3,msym) !vz_i
 89  real(dp),intent(inout) :: tnons(3,msym) !vz_i
 90 
 91 !Local variables ------------------------------
 92 ! intsym,inttn = intermediate real to swap the columns of the symmetry matrix
 93 ! bckbrvltt = backup to brvltt to compare the assigned and the input values
 94 !real(dp) :: tsec(2)
 95 !scalars
 96  integer :: bckbrvltt,ii,inversion,jj,kk,ierr
 97  integer :: intsym
 98  real(dp) :: inttn
 99  character(len=500) :: message
100 
101 ! *************************************************************************
102 
103 !List of the input parameters
104 !DEBUG
105 !write(std_out,*)' gensymspgr : enter with:'
106 !write(std_out,*)' spgroup = ',spgroup
107 !write(std_out,*)' spgaxor = ',spgaxor
108 !write(std_out,*)' spgorig = ',spgorig
109 !write(std_out,*)' brvltt  = ',brvltt
110 !ENDDEBUG
111 
112 !Assume that the value of spgroupma is consistent with the one of spgroup
113 !(this has been checked earlier)
114 
115 !Tests for consistency first the space group number and then the orientation
116 !Checks the space group number
117  if (.not.(spgroup>0 .and. spgroup<231) ) then
118    write(message, '(a,i0,a,a,a,a)' )&
119 &   'spgroup must be between 1 to 230, but is ',spgroup,ch10,&
120 &   'This is not allowed.  ',ch10,&
121 &   'Action: modify spgroup in the input file.'
122    ABI_ERROR(message)
123  end if
124 
125 !Checks the orientation
126  if (.not.(spgaxor>0 .and. spgaxor<10)) then
127    write(message, '(a,i12,a,a,a,a)' )&
128 &   'spgaxor must be from 1 to 9, but is',spgaxor,ch10,&
129 &   'This is not allowed.  ',ch10,&
130 &   'Action: modify spgaxor in the input file.'
131    ABI_ERROR(message)
132  end if
133 
134 !Checks the consistency between the origin and space group
135  if (spgorig==1 .or. spgorig==2) then
136  else
137    write(message, '(a,i4,a,a,a,a,a,a)' )&
138 &   'spgorig is',spgorig,ch10,&
139 &   'while it should be 0 or 1',ch10,&
140 &   'This is not allowed.  ',ch10,&
141 &   'Action: modify spgorig in the input file.'
142    ABI_ERROR(message)
143  end if
144 
145  if (spgorig>1) then
146    select case (spgroup)
147    case (48,50,59,68,70,85,86,88,125,126,129,130,133,134,137,138,&
148 &     141,142,201,203,222,224,227,228)
149    case default
150      write(message, '(a,a,a,a,a)' )&
151 &     'spgroup does not accept several origin choices',ch10,&
152 &     'This is not allowed.  ',ch10,&
153 &     'Action: modify spgorig in the input file.'
154      ABI_ERROR(message)
155    end select
156  end if
157 
158 !Checks for consistency between the orientation and space group
159  if (spgaxor>1) then
160    select case (spgroup)
161    case (3:74,146,148,155,160,161,166,167)
162    case default
163      write(message, '(a,a,a,a,a)' )&
164 &     'spgroup does not accept several orientations',ch10,&
165 &     'This is not allowed.  ',ch10,&
166 &     'Action: modify spgaxor or spgroup in the input file.'
167      ABI_ERROR(message)
168    end select
169  end if
170 
171  if (brvltt<-1 .or. brvltt>7)then
172    write(message, '(a,i4,a,a,a,a,a,a)' )&
173 &   'The input brvltt was ',brvltt,ch10,&
174 &   'and it should be an integer from -1 to 7',ch10,&
175 &   'This is not allowed.  ',ch10,&
176 &   'Action: modify brvltt  in the input file.'
177    ABI_ERROR(message)
178  end if
179 
180 !Assign nsym for each group according first to the order of the group
181 !Note that this value might be modified later:
182 !first because of the product with the inversion,
183 !second because of the centering operations
184  select case (spgroup)
185  case (1,2)
186    nsym=1
187  case (3:9)
188    nsym=2
189  case (143:148)
190    nsym=3
191  case (10:42,44:47,49,51:58,60:67,69,71:84,87)
192    nsym=4
193  case (149:176)
194    nsym=6
195  case (48,50,59,68,70,85,86,88:121,123,124,127,128,131,132,135,136,139,140)
196    nsym=8
197  case (177:200,202,204:206)
198    nsym=12
199  case (43,122,125,126,129,130,133,134,137,138)
200    nsym=16
201  case (201,207:219,221,223,225,226,229,230)
202    nsym=24
203  case (141,142)
204    nsym=32
205  case (203,220,222,224)
206    nsym=48
207  case (227,228)
208    nsym=192
209  end select
210 
211 !DEBUG
212 !write(std_out,*)'gensymspgr :  assigns nsym = ',nsym
213 !ENDDEBUG
214 
215 !Makes a backup to the brvltt for further comparison with the assigned value
216  bckbrvltt=brvltt
217 !Default brvltt
218  brvltt=1
219 
220 !call timab(47,1,tsec)
221 
222 !Assigns the first part of the symmetry operations:
223 !Rotation axis and mirror planes with or without translations,
224 !and sometimes also inversion operations.
225  select case (spgroup)
226  case (1:2)
227    symrel(:,:,1)=0
228    symrel(1,1,1)=1 ; symrel(2,2,1)=1 ; symrel(3,3,1)=1
229    tnons(:,1)=zero ; symafm(1)=1
230  case (3:15)
231    call symsgmono(brvltt,msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
232 &   spgroupma,symafm,symrel,tnons)
233  case (16:74)
234    call symsgortho(msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
235 &   spgroupma,symafm,symrel,tnons)
236  case (75:142)
237    call symsgtetra(msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
238 &   spgroupma,symafm,symrel,tnons)
239  case (143:194)
240    call symsghexa(brvltt,msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
241 &   spgroupma,symafm,symrel,tnons)
242  case (195:230)
243    call symsgcube(msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
244 &   spgroupma,symafm,symrel,tnons)
245  end select
246 
247 !call timab(47,2,tsec)
248 
249 !Assign the inversion center (if necessary).
250 !Note that for monoclinic space groups, the inversion was already
251 !assigned in symsgmono.f. Some other inversions have also been assigned in the
252 !corresponding system routine
253  inversion=0
254  select case (spgroup)
255  case (2,47,49,51:58,60:67,69,71:74,83,84,87,123,124,127,128,131,132,&
256 &   135,136,139,140,147,148,162:167,175,176,191:194,200,202,204:206,&
257 &   221,223,225,226,229,230)
258    inversion=1
259 !    Treat the magnetic part
260    if(shubnikov==3)then
261      select case (spgroup)
262      case(2) ! Triclinic
263        inversion=-1
264      case(47,49,51:58,60:67,69,71:74) ! Orthorhombic
265        select case (spgroupma)
266        case(251,253,259,261,267,268,271,279,280,283,291,292,293,297,307,&
267 &         308,309,313,323,324,325,329,339,340,341,345,355,356,359,367,368,371,&
268 &         379,380,381,385,395,396,399,407,408,411,419,420,421,425,435,437,443,&
269 &         444,445,449,459,460,461,465,471,472,473,477,483,484,487,493,494,497,&
270 &         503,504,507,513,514,517,523,525,529,531,535,537,541,542,545,550,552,556,557,560)
271          inversion=-1
272        end select
273      case(83,84,87,123,124,127,128,131,132,135,136,139,140) ! Tetragonal
274        select case (spgroupma)
275        case(46,47,54,55,62,63,70,71,78,79,84,85,341,344,346,347,353,356,358,359,365,&
276 &         368,370,371,377,380,382,383,389,392,394,395,401,404,406,407,&
277 &         413,416,418,419,425,428,430,431,437,440,442,443,449,452,454,455,&
278 &         461,464,466,467,473,476,478,479,485,488,490,491,497,500,502,503,&
279 &         509,512,514,515,521,524,526,527,533,536,538,539,543,546,548,549)
280          inversion=-1
281        end select
282      case(147,148,162:167,175,176,191:194) ! Hexagonal or rhombohedral
283        select case (spgroupma)
284        case(15,19,75,76,81,82,87,88,93,94,99,100,105,106,139,140,145,146,&
285 &         235,236,237,241,245,246,247,251,255,256,257,261,265,266,267,271)
286          inversion=-1
287        end select
288      case(200,202,204:206,221,223,225,226,229,230) ! Cubic
289        select case (spgroupma)
290        case(16,20,24,28,32,35,39,94,96,100,102,106,108,112,114,118,120,124,126,&
291 &         130,132,136,138,142,144,147,149)
292          inversion=-1
293        end select
294      end select
295    end if
296  end select
297 
298 !DEBUG
299 !write(std_out,*)' gensymspgr : before inversion'
300 !write(std_out,*)' Describe the different symmetry operations (index,symrel,tnons,symafm)'
301 !do ii=1,nsym
302 !write(std_out,'(i3,2x,9i3,3es12.2,i3)')ii,symrel(:,:,ii),tnons(:,ii),symafm(ii)
303 !end do
304 !ENDDEBUG
305 
306  if(inversion/=0)then
307    do ii=1,nsym        ! visit all the symmetries assigned before
308      do jj=1,3        ! visit the 3x3 matrix corresponding to the symmetry i
309        tnons(jj,nsym+ii)=-tnons(jj,ii)
310        do kk=1,3
311          symrel(jj,kk,nsym+ii)=-symrel(jj,kk,ii)
312        end do
313      end do
314      symafm(nsym+ii)=inversion*symafm(ii)
315    end do
316    nsym=nsym*2
317  end if
318 
319 !DEBUG
320 !write(std_out,*)' gensymspgr : after inversion'
321 !write(std_out,*)' Describe the different symmetry operations (index,symrel,tnons,symafm)'
322 !do ii=1,nsym
323 !write(std_out,'(i3,2x,9i3,3es12.2,i3)')ii,symrel(:,:,ii),tnons(:,ii),symafm(ii)
324 !end do
325 !ENDDEBUG
326 
327 !Assign the Bravais lattice to each space group to which it has not yet
328 !been assigned
329  select case (spgroup)
330  case (38:41)
331    brvltt=5                ! A
332  case (20,21,35:37,63:68)
333    brvltt=4                ! C
334  case (22,42,43,69,70,196,202,203,209,210,216,219,225:228)
335    brvltt=3                ! F
336  case (23,24,44:46,71:74,79,80,82,87,88,97,98,107:110,119:122,&
337 &   139:142,197,199,204,206,211,214,217,220,229,230)
338    brvltt=2                ! I
339  case (146,148,155,160,161,166,167)
340    if (spgaxor==1) then
341      brvltt=7
342    end if
343  end select
344 
345  if (bckbrvltt/=0 .and. bckbrvltt/=-1) then
346    if (bckbrvltt/=brvltt) then
347      write(message, '(a,i8,a,a,a,i8,a,a)' )&
348 &     'The assigned brvltt ',brvltt,' is not equal',ch10,&
349 &     'to the input value ',bckbrvltt,ch10,&
350 &     'Assume experienced user. Execution will continue.'
351      ABI_WARNING(message)
352    end if
353  end if
354 
355 !if(bckbrvltt>=0)then
356 !Complete the set of primitive symmetries by translations
357 !associated with brvltt.
358  select case (brvltt)
359 !  Bravais lattice type : A ! translation associated: b/2+c/2
360  case (5)
361    do ii=1,nsym
362      tnons(1,nsym+ii)=tnons(1,ii)
363      tnons(2,nsym+ii)=tnons(2,ii)+0.5
364      tnons(3,nsym+ii)=tnons(3,ii)+0.5
365      symrel(:,:,nsym+ii)=symrel(:,:,ii)
366      symafm(nsym+ii)=symafm(ii)
367    end do
368    nsym=nsym*2
369 
370 !    Bravais lattice type : B ! translation associated: a/2+c/2
371  case (6)
372    do ii=1,nsym
373      tnons(1,nsym+ii)=tnons(1,ii)+0.5
374      tnons(2,nsym+ii)=tnons(2,ii)
375      tnons(3,nsym+ii)=tnons(3,ii)+0.5
376      symrel(:,:,nsym+ii)=symrel(:,:,ii)
377      symafm(nsym+ii)=symafm(ii)
378    end do
379    nsym=nsym*2
380 
381 !    Bravais lattice type : C ! translation associated: a/2+b/2
382  case (4)
383    do ii=1,nsym
384      tnons(1,nsym+ii)=tnons(1,ii)+0.5
385      tnons(2,nsym+ii)=tnons(2,ii)+0.5
386      tnons(3,nsym+ii)=tnons(3,ii)
387      symrel(:,:,nsym+ii)=symrel(:,:,ii)
388      symafm(nsym+ii)=symafm(ii)
389    end do
390    nsym=nsym*2
391 
392 !    Bravais lattice type : F ! translations associated: a/2+b/2,b/2+c/2,b/2+c/2
393  case (3)
394 !    For space groups containing d elements, all the symmetry operations
395 !    have already been obtained
396     if(spgroup/=43 .and. spgroup/=227 .and. spgroup/=228)then  
397      do ii=1,nsym
398 !        First translation: a/2+b/2
399        tnons(1,nsym+ii)=tnons(1,ii)+0.5
400        tnons(2,nsym+ii)=tnons(2,ii)+0.5
401        tnons(3,nsym+ii)=tnons(3,ii)
402        symrel(:,:,nsym+ii)=symrel(:,:,ii)
403        symafm(nsym+ii)=symafm(ii)
404      end do
405 !      Second translation: b/2+c/2
406      do ii=1,nsym
407        tnons(1,2*nsym+ii)=tnons(1,ii)
408        tnons(2,2*nsym+ii)=tnons(2,ii)+0.5
409        tnons(3,2*nsym+ii)=tnons(3,ii)+0.5
410        symrel(:,:,2*nsym+ii)=symrel(:,:,ii)
411        symafm(2*nsym+ii)=symafm(ii)
412      end do
413 !      Third translation: a/2+c/2
414      do ii=1,nsym
415        tnons(1,3*nsym+ii)=tnons(1,ii)+0.5
416        tnons(2,3*nsym+ii)=tnons(2,ii)
417        tnons(3,3*nsym+ii)=tnons(3,ii)+0.5
418        symrel(:,:,3*nsym+ii)=symrel(:,:,ii)
419        symafm(3*nsym+ii)=symafm(ii)
420      end do
421      nsym=nsym*4
422    end if
423 
424 !    Bravais lattice type: I ! translation associated: a/2+b/2+c/2
425  case (2)
426 !    For space groups containing d elements, all the symmetry operations
427 !    have already been obtained
428    if(spgroup/=122 .and. spgroup/=141 .and. spgroup/=142 .and. &
429 &   spgroup/=220 )then
430      do ii=1,nsym        ! visit all the symmetries assigned before
431        tnons(:,nsym+ii)=tnons(:,ii)+0.5
432        symrel(:,:,nsym+ii)=symrel(:,:,ii)
433        symafm(nsym+ii)=symafm(ii)
434      end do
435      nsym=nsym*2
436    end if
437 
438 !    Bravais lattice type: R
439 !    translations for hexagonal axes ONLY: (2/3,1/3,1/3) & (1/3,2/3,2/3)
440 !    first translation (2/3,1/3,1/3)
441  case (7)
442    do ii=1,nsym
443      tnons(1,nsym+ii)=tnons(1,ii)+two_thirds
444      tnons(2,nsym+ii)=tnons(2,ii)+third
445      tnons(3,nsym+ii)=tnons(3,ii)+third
446      symrel(:,:,nsym+ii)=symrel(:,:,ii)
447      symafm(nsym+ii)=symafm(ii)
448    end do
449 !    Second translation (1/3,2/3,2/3)
450    do ii=1,nsym
451      tnons(1,2*nsym+ii)=tnons(1,ii)+third
452      tnons(2,2*nsym+ii)=tnons(2,ii)+two_thirds
453      tnons(3,2*nsym+ii)=tnons(3,ii)+two_thirds
454      symrel(:,:,2*nsym+ii)=symrel(:,:,ii)
455      symafm(2*nsym+ii)=symafm(ii)
456    end do
457    nsym=nsym*3
458 
459  end select
460 
461 !end if
462 
463 !Translate tnons in the ]-0.5,0.5] interval
464  tnons(:,1:nsym)=tnons(:,1:nsym)-nint(tnons(:,1:nsym)-1.0d-8)
465 
466 !Orientations for the orthorhombic space groups
467 !WARNING : XG 000620 : I am not sure that this coding is correct !!
468  if (spgroup>15 .and. spgroup <75) then
469    select case (spgaxor)
470    case (1)             ! abc
471      write(std_out,*)' the choosen orientation corresponds to: abc; the proper one'
472    case (2)             ! cab
473      do ii=1,nsym
474        intsym=symrel(1,1,ii)
475        symrel(1,1,ii)=symrel(2,2,ii)
476        symrel(2,2,ii)=intsym
477        inttn=tnons(1,ii)
478        tnons(1,ii)=tnons(2,ii)
479        tnons(2,ii)=inttn
480      end do
481      do ii=1,nsym
482        intsym=symrel(1,1,ii)
483        symrel(1,1,ii)=symrel(3,3,ii)
484        symrel(3,3,ii)=intsym
485        inttn=tnons(1,ii)
486        tnons(1,ii)=tnons(3,ii)
487        tnons(3,ii)=inttn
488      end do
489      write(std_out,*)' the choosen orientation corresponds to:  cab'
490    case (3)             ! bca
491      do ii=1,nsym
492        intsym=symrel(1,1,ii)
493        symrel(1,1,ii)=symrel(2,2,ii)
494        symrel(2,2,ii)=intsym
495        inttn=tnons(1,ii)
496        tnons(1,ii)=tnons(2,ii)
497        tnons(2,ii)=inttn
498      end do
499      do ii=1,nsym
500        intsym=symrel(2,2,ii)
501        symrel(2,2,ii)=symrel(3,3,ii)
502        symrel(3,3,ii)=intsym
503        inttn=tnons(2,ii)
504        tnons(2,ii)=tnons(3,ii)
505        tnons(3,ii)=inttn
506      end do
507      write(std_out,*)' the choosen orientation corresponds to:  bca'
508    case (4)             ! acb
509      do ii=1,nsym
510        intsym=symrel(2,2,ii)
511        symrel(2,2,ii)=symrel(3,3,ii)
512        symrel(3,3,ii)=intsym
513        inttn=tnons(1,ii)
514        tnons(2,ii)=tnons(3,ii)
515        tnons(3,ii)=inttn
516      end do
517      write(std_out,*)' the choosen orientation corresponds to:  acb'
518    case (5)             ! bac
519      do ii=1,nsym
520        intsym=symrel(1,1,ii)
521        symrel(1,1,ii)=symrel(2,2,ii)
522        symrel(2,2,ii)=intsym
523        inttn=tnons(1,ii)
524        tnons(1,ii)=tnons(2,ii)
525        tnons(2,ii)=inttn
526      end do
527      write(std_out,*)' the choosen orientation corresponds to:  bac'
528    case (6)             ! cba
529      do ii=1,nsym
530        intsym=symrel(1,1,ii)
531        symrel(1,1,ii)=symrel(3,3,ii)
532        symrel(3,3,ii)=intsym
533        inttn=tnons(1,ii)
534        tnons(1,ii)=tnons(3,ii)
535        tnons(3,ii)=inttn
536      end do
537      write(std_out,*)' the choosen orientation corresponds to:  cba'
538    end select
539  end if
540 
541 !DEBUG
542 !write(std_out,*)' gensymspgr  : out of the Bravais lattice, nsym is',nsym
543 !ENDDEBUG
544 
545  call sg_multable(nsym,symafm,symrel,ierr,tnons=tnons,tnons_tol=tol5)
546  if (ierr/=0) then
547    call print_symmetries(nsym,symrel,tnons,symafm)
548  end if
549 
550  ABI_CHECK(ierr==0,"Error in group closure")
551 
552 end subroutine gensymspgr