TABLE OF CONTENTS


ABINIT/prtspgroup [ Functions ]

[ Top ] [ Functions ]

NAME

 prtspgroup

FUNCTION

 Print the space group (first, the dataset)

COPYRIGHT

 Copyright (C) 2002-2018 ABINIT group (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 .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt.

INPUTS

  bravais(11)=characteristics of Bravais lattice (see symlatt.f)
  genafm(3)=generator of magnetic translations, in case of
            Shubnikov type IV magnetic groups (if zero, the group is
            not a type IV magnetic group)
  iout=unit number of output file
  jdtset= actual number of the dataset to be read
  ptgroupma=magnetic point group, in case of
            Shubnikov type III magnetic groups (if zero, the group is
            not a type III magnetic group)
  spgroup=space group number

OUTPUT

PARENTS

      memory_eval

CHILDREN

      ptgmadata,spgdata,wrtout,xred2xcart

SOURCE

 38 #if defined HAVE_CONFIG_H
 39 #include "config.h"
 40 #endif
 41 
 42 #include "abi_common.h"
 43 
 44 
 45 subroutine prtspgroup(bravais,genafm,iout,jdtset,ptgroupma,spgroup)
 46 
 47  use defs_basis
 48  use m_errors
 49  use m_profiling_abi
 50 
 51  use m_geometry,  only : xred2xcart
 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 'prtspgroup'
 57  use interfaces_14_hidewrite
 58  use interfaces_41_geometry, except_this_one => prtspgroup
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  integer,intent(in) :: iout,jdtset,ptgroupma,spgroup
 66 !arrays
 67  integer,intent(in) :: bravais(11)
 68  real(dp),intent(inout) :: genafm(3)
 69 
 70 !Local variables -------------------------------
 71 !scalars
 72  integer :: center,iholohedry,ii,shubnikov,spgaxor,spgorig,sporder,sumgen
 73  character(len=1) :: brvsb
 74  character(len=10) :: ptgrpmasb
 75  character(len=15) :: intsb,ptintsb,ptschsb,schsb
 76  character(len=35) :: intsbl
 77  character(len=500) :: message
 78  character(len=80) :: bravais_name
 79 !arrays
 80  integer :: genafmint(3)
 81  real(dp) :: genafmconv(3),rprimdconv(3,3)
 82 
 83 !*************************************************************************
 84 
 85 !DEBUG
 86 !write(std_out,*)' prtspgroup : enter '
 87 !write(std_out,*)' ptgroupma=',ptgroupma
 88 !write(std_out,*)' genafm(:)=',genafm(:)
 89 !ENDDEBUG
 90 
 91  center=bravais(2)
 92  iholohedry=bravais(1)
 93 
 94 !Determine the magnetic type
 95  shubnikov=1
 96  if(ptgroupma/=0)shubnikov=3
 97  if(sum(abs(genafm(:)))>tol6)then
 98    shubnikov=4
 99 !  Produce genafm in conventional axes,
100    rprimdconv(:,1)=bravais(3:5)
101    rprimdconv(:,2)=bravais(6:8)
102    rprimdconv(:,3)=bravais(9:11)
103    if(center/=0)rprimdconv(:,:)=rprimdconv(:,:)*half
104    call xred2xcart(1,rprimdconv,genafmconv,genafm)
105 !  Gives the associated translation, with components in the
106 !  interval ]-0.5,0.5] .
107    genafmconv(:)=genafmconv(:)-nint(genafmconv(:)-tol6)
108    do ii=1,3
109      genafmint(ii)=-1
110      if(abs(genafmconv(ii)-zero)<tol6)genafmint(ii)=0
111      if(abs(genafmconv(ii)-half)<tol6)genafmint(ii)=1
112    end do
113    if(minval(genafmint(:))==-1)then
114      write(message, '(3a,3es12.2,a)' )&
115 &     'The magnetic translation generator,',ch10,&
116 &     'genafmconv(:)=',genafmconv(:),&
117 &     'could not be identified.'
118      MSG_BUG(message)
119    end if
120  end if
121 
122 !Determine whether the space group can be printed
123  if(iholohedry<=0)then
124    if(jdtset/=0)then
125      write(message,'(a,a,i5,a)')ch10,&
126 &     ' DATASET',jdtset,' : the unit cell is not primitive'
127    else
128      write(message,'(a,a)')ch10,&
129 &     ' Symmetries : the unit cell is not primitive'
130    end if
131    call wrtout(std_out,message,'COLL')
132    call wrtout(iout,message,'COLL')
133  else if(spgroup==0)then
134    if(jdtset/=0)then
135      write(message,'(a,a,i5,a)')ch10,&
136 &     ' DATASET',jdtset,' : the space group has not been recognized'
137    else
138      write(message,'(a,a)')ch10,&
139 &     ' Symmetries : the space group has not been recognized'
140    end if
141    call wrtout(std_out,message,'COLL')
142    call wrtout(iout,message,'COLL')
143  else
144 
145 !  ------------------------------------------------------------------
146 !  The space group can be printed
147 
148 !  Determine the Bravais lattice
149 
150    bravais_name=' (the Bravais lattice could not be identified)'
151 
152    if(iholohedry==7)then ! Cubic
153 
154      if(center==0) then
155        if(shubnikov/=4)bravais_name='cP (primitive cubic)'
156        if(shubnikov==4)bravais_name='cP_I (primitive cubic, inner magnetic, #33)'
157      else if(center==-1) then
158        if(shubnikov/=4)bravais_name='cI (body-center cubic)' ! Only non-magnetic is possible
159      else if(center==-3) then
160        if(shubnikov/=4)bravais_name='cF (face-center cubic)'
161        if(shubnikov==4)bravais_name='cF_s (face-center cubic, simple cubic magnetic, #35)'
162      end if
163 
164    else if(iholohedry==4)then ! Tetragonal
165 
166      if(center==0) then
167        if(shubnikov/=4)bravais_name='tP (primitive tetrag.)'
168        if(shubnikov==4)then
169          sumgen=sum(genafmint(:))
170          if(sumgen==1)bravais_name='tP_c (primitive tetrag., c-magnetic, #23)'
171          if(sumgen==2)bravais_name='tP_C (primitive tetrag., C-magnetic, #24)'
172          if(sumgen==3)bravais_name='tP_I (primitive tetrag., centered magnetic, #25)'
173        end if
174      else if(center==-1)then
175        if(shubnikov/=4)bravais_name='tI (body-center tetrag.)'
176        if(shubnikov==4)bravais_name='tI_c (body-center tetrag., simple tetragonal magnetic, #27)'
177      end if
178 
179    else if(iholohedry==3)then ! Orthorhombic
180 
181      if(center==0) then
182        if(shubnikov/=4)bravais_name='oP (primitive ortho.)'
183        if(shubnikov==4)then
184          sumgen=sum(genafmint(:))
185          if(sumgen==1)then
186            if(genafmint(1)==1)bravais_name='oP_a (primitive ortho., a-magnetic, #11)'
187            if(genafmint(2)==1)bravais_name='oP_b (primitive ortho., b-magnetic, #11)'
188            if(genafmint(3)==1)bravais_name='oP_c (primitive ortho., c-magnetic, #11)'
189          else if(sumgen==2)then
190            if(genafmint(1)==0)bravais_name='oP_A (primitive ortho., A-magnetic, #12)'
191            if(genafmint(2)==0)bravais_name='oP_B (primitive ortho., B-magnetic, #12)'
192            if(genafmint(3)==0)bravais_name='oP_C (primitive ortho., C-magnetic, #12)'
193          else if(sumgen==3)then
194            bravais_name='oP_I (primitive ortho., centered magnetic, #13)'
195          end if
196        end if
197      else if(center==-1)then
198        if(shubnikov/=4)bravais_name='oI (body-center ortho.)'
199        if(shubnikov==4)bravais_name='oI_c (body-center ortho., simple ortho. magn., #21)'
200      else if(center==1 .or. center==2 .or. center==3)then
201        if(shubnikov/=4) bravais_name='oC (1-face-center ortho.)'
202        if(shubnikov==4)then
203          sumgen=sum(genafmint(:))
204          if(sumgen==1)then
205 !          One should work more to distinguish these magnetic groups
206            bravais_name='oC_(a,b,c) (1-face-cent. ortho., 1-magn., #15 or 16)'
207          else if(sumgen==2)then
208            bravais_name='oC_A (1-face-centered ortho., 1-face-magnetic, #17)'
209          else if(sumgen==3)then
210            bravais_name='oC_c (C-face-centered ortho., c-magnetic, #15)'
211          end if
212        end if
213      else if(center==-3)then
214        if(shubnikov/=4)bravais_name='oF (face-center ortho.)'
215        if(shubnikov==4)bravais_name='oF_s (face-center ortho., simple ortho. magnetic, #19)'
216      end if
217 
218    else if(iholohedry==6)then ! Hexagonal
219 
220      if(shubnikov/=4)bravais_name='hP (primitive hexag.)'
221      if(shubnikov==4)bravais_name='hP_c (primitive hexag., c-magnetic, #29)'
222 
223    else if(iholohedry==5)then ! Rhombohedral
224 
225      if(shubnikov/=4)bravais_name='hR (rhombohedral)'
226      if(shubnikov==4)bravais_name='hR_I (rhombohedral, centered magnetic, #31)'
227 
228    else if(iholohedry==2)then ! Monoclinic
229 
230      if(center==0)then
231        if(shubnikov/=4)bravais_name='mP (primitive monocl.)'
232        if(shubnikov==4)then
233          sumgen=sum(genafmint(:))
234          if(sumgen==1)then
235            if(genafmint(1)==1)bravais_name='mP_a (primitive monocl., a-magnetic, #5)'
236            if(genafmint(2)==1)bravais_name='mP_b (primitive monocl., b-magnetic, #4)'
237            if(genafmint(3)==1)bravais_name='mP_c (primitive monocl., c-magnetic, #5)'
238          else if(sumgen==2)then
239            if(genafmint(1)==0)bravais_name='mP_A (primitive monocl., A-magnetic, #6)'
240            if(genafmint(2)==0)bravais_name='mP_B (primitive monocl., B-magnetic, #6)'
241            if(genafmint(3)==0)bravais_name='mP_C (primitive monocl., C-magnetic, #6)'
242          end if
243        end if
244      else if(center==3)then
245        if(shubnikov/=4)bravais_name='mC (1-face-center monocl.)'
246        if(shubnikov==4)then
247          if(genafmint(3)==1)bravais_name='mC_c (C-face-center monocl., c-magnetic, #8)'
248          if(genafmint(3)/=1)bravais_name='mC_a (C-face-center monocl., a-magnetic, #9)'
249        end if
250      else if(center==-3)then
251        if(shubnikov/=4)bravais_name='(reduction of face-center)'
252      end if
253 
254    else if(iholohedry==1)then ! Triclinic
255 
256      if(shubnikov/=4)bravais_name='aP (primitive triclinic)'
257      if(shubnikov==4)bravais_name='aP_s (primitive triclinic, simple magnetic, #2)'
258 
259    end if
260 
261 !  Determine the symbol of the Fedorov space group
262    spgaxor=1 ; spgorig=1
263    call spgdata(brvsb,intsb,intsbl,ptintsb,ptschsb,&
264 &   schsb,spgaxor,spgroup,sporder,spgorig)
265 
266 !  Prepare print of the dataset, symmetry point group, Bravais lattice
267    if(shubnikov==1)then
268 
269      if(jdtset/=0)then
270        write(message,'(a,a,i5,a,a,a,a,i3,a,a,a)' )ch10,&
271 &       ' DATASET',jdtset,' : space group ',trim(brvsb),trim(intsb),' (#',spgroup,')',&
272 &       '; Bravais ',trim(bravais_name)
273      else
274        write(message,'(a,a,a,a,a,i3,a,a,a)' )ch10,&
275 &       ' Symmetries : space group ',trim(brvsb),trim(intsb),' (#',spgroup,')',&
276 &       '; Bravais ',trim(bravais_name)
277      end if
278      call wrtout(std_out,message,'COLL')
279      call wrtout(iout,message,'COLL')
280 
281    else if(shubnikov==3)then
282 
283      if(jdtset/=0)then
284        write(message,'(a,a,i5,a)' )ch10,&
285 &       ' DATASET',jdtset,' : magnetic group, Shubnikov type III '
286      else
287        write(message,'(2a)' )ch10,&
288 &       ' Magnetic group, Shubnikov type III '
289      end if
290      call wrtout(std_out,message,'COLL')
291      call wrtout(iout,message,'COLL')
292 
293      write(message,'(a,a,a,a,i3,a,a,a)' )&
294 &     ' Fedorov space group ',trim(brvsb),trim(intsb),' (#',spgroup,')',&
295 &     '; Bravais ',trim(bravais_name)
296      call wrtout(std_out,message,'COLL')
297      call wrtout(iout,message,'COLL')
298 
299      call ptgmadata(ptgroupma,ptgrpmasb)
300 
301      write(message,'(3a,i3,a)' )&
302 &     ' Magnetic point group ',trim(ptgrpmasb),' (#',ptgroupma,')'
303      call wrtout(std_out,message,'COLL')
304      call wrtout(iout,message,'COLL')
305 
306    else if(shubnikov==4)then
307 
308      if(jdtset/=0)then
309        write(message,'(a,a,i5,a)' )ch10,&
310 &       ' DATASET',jdtset,' : magnetic group, Shubnikov type IV '
311      else
312        write(message,'(2a)' )ch10,&
313 &       ' Magnetic group, Shubnikov type IV '
314      end if
315      call wrtout(std_out,message,'COLL')
316      call wrtout(iout,message,'COLL')
317 
318      write(message,'(a,a,a,a,i3,a)' )&
319 &     ' Fedorov space group ',trim(brvsb),trim(intsb),' (#',spgroup,')'
320      call wrtout(std_out,message,'COLL')
321      call wrtout(iout,message,'COLL')
322 
323      write(message,'(2a)' )&
324 &     ' Magnetic Bravais lattice ',trim(bravais_name)
325      call wrtout(std_out,message,'COLL')
326      call wrtout(iout,message,'COLL')
327 
328    end if
329  end if
330 
331 end subroutine prtspgroup