1!{\src2tex{textfont=tt}}
2!!****f* ABINIT/symlist_bcc
3!! NAME
4!! symlist_bcc
5!!
6!! FUNCTION
7!! Determine the space group from the number and type of symmetry operations
8!! BCC case.
9!!
10!! COPYRIGHT
11!! Copyright (C) 2000-2016 ABINIT group (RC)
12!! This file is distributed under the terms of the
13!! GNU General Public License, see ~abinit/COPYING
14!! or http://www.gnu.org/copyleft/gpl.txt .
15!! For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
16!!
17!! INPUTS
18!! additional_info=information that is needed beyond n_axes, in order
19!!  to discriminate between specific space groups
20!! nsym=actual number of symmetries
21!! n_axes(31)=array containing the number of all the possible symmetry operations
22!!
23!! OUTPUT
24!! spgroup=space group number ; returns 0 if not found
25!!
26!! NOTES
27!!
28!! The list of symmetry operations is for the conventional cell
29!!
30!! TODO
31!! For the time being there are several groups where uncertainties still exist
32!! This will be solved in the very next ABINIT version
33!!
34!! PARENTS
35!!      symspgr
36!!
37!! CHILDREN
38!!
39!! SOURCE
40
41#if defined HAVE_CONFIG_H
42#include "config.h"
43#endif
44
45#include "abi_common.h"
46
47
48subroutine symlist_bcc(additional_info,nsym,n_axes,spgroup)
49
50 use defs_basis
51 use m_profiling_abi
52
53!This section has been created automatically by the script Abilint (TD).
54!Do not modify the following lines by hand.
55#undef ABI_FUNC
56#define ABI_FUNC 'symlist_bcc'
57!End of the abilint section
58
59 implicit none
60
61!Arguments ------------------------------------
62!scalars
63 integer,intent(in) :: additional_info,nsym
64 integer,intent(out) :: spgroup
65!arrays
66 integer,intent(in) :: n_axes(31)
67
68!Local variables-------------------------------
69!character(len=500) :: message
70!scalars
71 integer :: ii
72!arrays
73 integer :: n_axest(31)
74
75!**************************************************************************
76
77!DEBUG
78!write(std_out,*) ' symlist_bcc : enter '
79!write(std_out,*) ' nsym = ', nsym
80!write(std_out,'(a,10i3)' ) ' n_axes(1:10) =',n_axes(1:10)
81!write(std_out,'(a,10i3)' ) ' n_axes(11:20)=',n_axes(11:20)
82!write(std_out,'(a,11i3)' ) ' n_axes(21:31)=',n_axes(21:31)
83!ENDDEBUG
84
85 spgroup=0
86
87 select case(nsym)
88
89 case(4)
90!    XG021015, from LSi.
91!    The coding of this case is different because of a compiler problem on
92!    SUN. Seems that an internal limit is reached. Do not know why.
93   do ii=1,3
94     if(ii==1)then
95       n_axest(:)=0 ; n_axest(7)=1 ; n_axest(8)=1 ; n_axest(9)=1 ; n_axest(20)=1
96       spgroup=5
97     else if(ii==2)then
98       n_axest(:)=0 ; n_axest(7)=1 ; n_axest(8)=1 ; n_axest(15)=1 ; n_axest(18)=1
99       spgroup=8
100     else if(ii==3)then
101       n_axest(:)=0 ; n_axest(7)=1 ; n_axest(8)=1 ; n_axest(16)=1 ; n_axest(18)=1
102       spgroup=9
103     end if
104     if(sum((n_axes-n_axest)**2)==0)exit
105     spgroup=0
106   end do
107
108 case(8)
109
110   n_axest=(/0,0,0,0,2,0,1,1,1,0,  0,0,0,0,0,2,0,0,0,1,  0,0,0,0,0,0,0,0,0,0,0/)
111   if(sum((n_axes-n_axest)**2)==0) spgroup=15
112   n_axest=(/0,0,0,0,0,0,1,1,3,0,  0,0,0,0,0,0,0,0,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
113   if(sum((n_axes-n_axest)**2)==0) then
114     if(additional_info==1) spgroup=23
115     if(additional_info==2) spgroup=24
116   end if
117   n_axest=(/0,0,0,0,0,0,1,1,1,0,  0,0,0,0,2,0,0,2,0,1,  0,0,0,0,0,0,0,0,0,0,0/)
118   if(sum((n_axes-n_axest)**2)==0) spgroup=44
119   n_axest=(/0,0,0,0,0,0,1,1,1,0,  0,0,0,0,0,4,0,0,0,1,  0,0,0,0,0,0,0,0,0,0,0/)
120   if(sum((n_axes-n_axest)**2)==0) spgroup=45
121   n_axest=(/0,0,0,0,0,0,1,1,1,0,  0,0,0,0,1,2,0,1,0,1,  0,0,0,0,0,0,0,0,0,0,0/)
122   if(sum((n_axes-n_axest)**2)==0) spgroup=46
123   n_axest=(/0,0,0,0,0,0,1,1,2,0,  0,2,0,0,0,0,0,0,0,0,  0,0,0,0,2,0,0,0,0,0,0/)
124   if(sum((n_axes-n_axest)**2)==0) spgroup=79
125   n_axest=(/0,0,0,0,0,0,1,1,2,0,  0,0,0,0,0,0,0,0,0,0,  0,0,0,4,0,0,0,0,0,0,0/)
126   if(sum((n_axes-n_axest)**2)==0) spgroup=80
127
128   n_axest=(/0,4,0,0,0,0,1,1,2,0,  0,0,0,0,0,0,0,0,0,0,  0,0,0,0,0,0,0,0,0,0,0/)
129   if(sum((n_axes-n_axest)**2)==0) spgroup=82
130
131 case(16)
132
133   n_axest=(/0,0,0,0,2,0,1,1,3,0,  0,0,0,0,3,0,0,3,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
134   if(sum((n_axes-n_axest)**2)==0) spgroup=71
135   n_axest=(/0,0,0,0,2,0,1,1,3,0,  0,0,0,0,1,4,0,1,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
136   if(sum((n_axes-n_axest)**2)==0) spgroup=72
137   n_axest=(/0,0,0,0,2,0,1,1,3,0,  0,0,0,0,0,6,0,0,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
138   if(sum((n_axes-n_axest)**2)==0) spgroup=73
139   n_axest=(/0,0,0,0,2,0,1,1,3,0,  0,0,0,0,2,2,0,2,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
140   if(sum((n_axes-n_axest)**2)==0) spgroup=74
141   n_axest=(/0,0,0,0,0,0,1,1,6,0,  0,2,0,0,0,0,0,0,0,0,  4,0,0,0,2,0,0,0,0,0,0/)
142   if(sum((n_axes-n_axest)**2)==0) spgroup=97
143   n_axest=(/0,0,0,0,0,0,1,1,6,0,  0,0,0,0,0,0,0,0,0,0,  4,0,0,4,0,0,0,0,0,0,0/)
144   if(sum((n_axes-n_axest)**2)==0) spgroup=98
145   n_axest=(/0,0,0,0,0,0,1,1,2,0,  0,2,0,0,8,0,0,0,0,0,  0,0,0,0,2,0,0,0,0,0,0/)
146   if(sum((n_axes-n_axest)**2)==0) spgroup=107
147   n_axest=(/0,0,0,0,0,0,1,1,2,0,  0,2,0,0,4,4,0,0,0,0,  0,0,0,0,2,0,0,0,0,0,0/)
148   if(sum((n_axes-n_axest)**2)==0) spgroup=108
149   n_axest=(/0,0,0,0,0,0,1,1,2,0,  0,0,0,0,4,0,4,0,0,0,  0,0,0,4,0,0,0,0,0,0,0/)
150   if(sum((n_axes-n_axest)**2)==0) spgroup=109
151   n_axest=(/0,0,0,0,0,0,1,1,2,0,  0,0,0,0,0,4,4,0,0,0,  0,0,0,4,0,0,0,0,0,0,0/)
152   if(sum((n_axes-n_axest)**2)==0) spgroup=110
153
154   n_axest=(/0,4,0,0,2,0,1,1,2,0,  0,2,0,0,2,0,0,0,0,0,  0,0,0,0,2,0,0,0,0,0,0/)
155   if(sum((n_axes-n_axest)**2)==0) spgroup=87
156   n_axest=(/0,4,0,0,2,0,1,1,2,0,  0,0,0,0,0,2,0,0,0,0,  0,0,0,4,0,0,0,0,0,0,0/)
157   if(sum((n_axes-n_axest)**2)==0) spgroup=88
158   n_axest=(/0,4,0,0,0,0,1,1,2,0,  0,0,0,0,4,0,0,0,0,0,  4,0,0,0,0,0,0,0,0,0,0/)
159   if(sum((n_axes-n_axest)**2)==0) spgroup=119
160   n_axest=(/0,4,0,0,0,0,1,1,2,0,  0,0,0,0,0,4,0,0,0,0,  4,0,0,0,0,0,0,0,0,0,0/)
161   if(sum((n_axes-n_axest)**2)==0) spgroup=120
162   n_axest=(/0,4,0,0,0,0,1,1,6,0,  0,0,0,0,4,0,0,0,0,0,  0,0,0,0,0,0,0,0,0,0,0/)
163   if(sum((n_axes-n_axest)**2)==0) spgroup=121
164   n_axest=(/0,4,0,0,0,0,1,1,6,0,  0,0,0,0,0,0,4,0,0,0,  0,0,0,0,0,0,0,0,0,0,0/)
165   if(sum((n_axes-n_axest)**2)==0) spgroup=122
166
167 case(24)
168
169   n_axest=(/0,0,0,0,0,0,1,1,3,16, 0,0,0,0,0,0,0,0,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
170   if(sum((n_axes-n_axest)**2)==0) then
171     if(additional_info==1) spgroup=197
172     if(additional_info==2) spgroup=199
173   end if
174
175 case(32)
176
177   n_axest=(/0,4,0,0,2,0,1,1,6,0,  0,2,0,0,10,0,0,0,0,0, 4,0,0,0,2,0,0,0,0,0,0/)
178   if(sum((n_axes-n_axest)**2)==0) spgroup=139
179   n_axest=(/0,4,0,0,2,0,1,1,6,0,  0,2,0,0,6,4,0,0,0,0,  4,0,0,0,2,0,0,0,0,0,0/)
180   if(sum((n_axes-n_axest)**2)==0) spgroup=140
181   n_axest=(/0,4,0,0,2,0,1,1,6,0,  0,0,0,0,4,2,4,0,0,0,  4,0,0,4,0,0,0,0,0,0,0/)
182   if(sum((n_axes-n_axest)**2)==0) spgroup=141
183   n_axest=(/0,4,0,0,2,0,1,1,6,0,  0,0,0,0,0,6,4,0,0,0,  4,0,0,4,0,0,0,0,0,0,0/)
184   if(sum((n_axes-n_axest)**2)==0) spgroup=142
185
186 case(48)
187
188   n_axest=(/0,0,16,0,2,0,1,1,3,16,  0,0,0,0,3,0,0,3,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
189   if(sum((n_axes-n_axest)**2)==0) spgroup=204
190   n_axest=(/0,0,16,0,2,0,1,1,3,16,  0,0,0,0,0,0,0,6,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
191   if(sum((n_axes-n_axest)**2)==0) spgroup=206
192   n_axest=(/0,0,0,0,0,0,1,1,12,16,  0,6,0,0,0,0,0,0,0,6,  0,0,0,0,6,0,0,0,0,0,0/)
193   if(sum((n_axes-n_axest)**2)==0) spgroup=211
194   n_axest=(/0,0,0,0,0,0,1,1,12,16,  0,0,0,0,0,0,0,0,0,6,  0,0,0,12,0,0,0,0,0,0,0/)
195   if(sum((n_axes-n_axest)**2)==0) spgroup=214
196   n_axest=(/0,12,0,0,0,0,1,1,3,16,  0,0,0,0,6,0,0,6,0,3,  0,0,0,0,0,0,0,0,0,0,0/)
197   if(sum((n_axes-n_axest)**2)==0) spgroup=217
198   n_axest=(/0,12,0,0,0,0,1,1,3,16,  0,0,0,0,0,0,12,0,0,3, 0,0,0,0,0,0,0,0,0,0,0/)
199   if(sum((n_axes-n_axest)**2)==0) spgroup=220
200
201 case(96)
202
203!    Note that the identification of the mirror planes is still ambiguous for cI
204   n_axest=(/0,12,16,0,2,0,1,1,12,16,  0,6,0,0,9,0,0,9,0,6,  0,0,0,0,6,0,0,0,0,0,0/)
205   if(sum((n_axes-n_axest)**2)==0) spgroup=229
206   n_axest=(/0,12,16,0,2,0,1,1,12,16,  0,0,0,0,0,0,12,6,0,6, 0,0,0,12,0,0,0,0,0,0,0/)
207   if(sum((n_axes-n_axest)**2)==0) spgroup=230
208
209 end select
210
211end subroutine symlist_bcc
212!!***
213