TABLE OF CONTENTS


ABINIT/m_ab7_symmetry [ Modules ]

[ Top ] [ Modules ]

NAME

 m_ab7_symmetry

FUNCTION

COPYRIGHT

  Copyright (C) 2008-2024 ABINIT group (DC)
  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

  15 #if defined HAVE_CONFIG_H
  16 #include "config.h"
  17 #endif
  18 
  19 #include "abi_common.h"
  20 
  21 module m_ab7_symmetry
  22 
  23   use defs_basis
  24   use m_abicore
  25 
  26   use m_symtk,     only : mati3inv, mati3det, symatm, symcharac
  27   use m_symfind,   only : symfind, symanal, symlatt
  28   use m_geometry,  only : metric
  29   use m_spgdata,   only : spgdata
  30 
  31   implicit none
  32 
  33   private
  34 
  35   integer, parameter, public :: AB7_MAX_SYMMETRIES = 384
  36 
  37   type, public :: symmetry_type
  38      ! The input characteristics
  39      real(dp) :: tolsym
  40      real(dp) :: rprimd(3,3), gprimd(3,3), rmet(3,3)
  41      integer :: nAtoms
  42      integer, pointer :: typeAt(:)
  43      real(dp), pointer :: xRed(:,:)
  44 
  45      logical :: withField
  46      real(dp) :: field(3)
  47 
  48      logical :: withJellium
  49 
  50      integer :: nzchempot
  51 
  52      integer :: withSpin
  53      real(dp), pointer :: spinAt(:,:)
  54 
  55      logical :: withSpinOrbit
  56 
  57      integer :: vacuum(3)
  58 
  59      ! The output characteristics
  60      ! The bravais parameters
  61      integer :: nBravSym
  62      integer :: bravais(11), bravSym(3, 3, AB7_MAX_SYMMETRIES)
  63      ! The symmetry matrices
  64      logical  :: auto
  65      integer  :: nSym
  66      integer, pointer  :: sym(:,:,:)
  67      real(dp), pointer :: transNon(:,:)
  68      integer, pointer  :: symAfm(:)
  69      ! Some additional information
  70      integer          :: multiplicity
  71      real(dp)         :: genAfm(3)
  72      integer          :: spaceGroup, pointGroupMagn
  73      integer, pointer :: indexingAtoms(:,:,:)
  74   end type symmetry_type
  75 
  76   ! We store here a list of symmetry objects to be able to
  77   ! call several symmetry operations on different objects.
  78   ! The simplest portable way to do it, is to create
  79   ! a list of Fortran structure and to use the list index
  80   ! as an identifier that can be given to the other languages.
  81   type, private :: symmetry_list
  82      integer                       :: id
  83      type(symmetry_list),  pointer :: next
  84      type(symmetry_type)           :: data
  85   end type symmetry_list
  86   type(symmetry_list), pointer :: my_symmetries
  87   integer :: n_symmetries = 0
  88 
  89   logical, private, parameter :: AB_DBG = .false.
  90 
  91   public :: symmetry_new
  92   public :: symmetry_free
  93   public :: symmetry_set_tolerance
  94   public :: symmetry_set_lattice
  95   public :: symmetry_set_structure
  96   public :: symmetry_set_collinear_spin
  97   public :: symmetry_set_spin
  98   public :: symmetry_set_spin_orbit
  99   public :: symmetry_set_field
 100   public :: symmetry_set_jellium
 101   public :: symmetry_set_periodicity
 102   public :: symmetry_set_n_sym
 103 
 104   public :: symmetry_get_from_id
 105   public :: symmetry_get_n_atoms
 106   public :: symmetry_get_n_sym
 107   public :: symmetry_get_multiplicity
 108   public :: symmetry_get_bravais
 109   public :: symmetry_get_matrices
 110   public :: symmetry_get_matrices_p
 111   public :: symmetry_get_group
 112   public :: symmetry_get_equivalent_atom
 113   public :: symmetry_get_type
 114 
 115 contains
 116 
 117   subroutine new_item(token)
 118 
 119     type(symmetry_list), pointer :: token
 120 
 121     ! We allocate a new list token and prepend it.
 122     if (AB_DBG) write(std_err,*) "AB symmetry: create a new token."
 123 
 124     ! Init case, very first call.
 125     if (n_symmetries == 0) then
 126        nullify(my_symmetries)
 127     end if
 128 
 129     ! Normal treatment.
 130     n_symmetries = n_symmetries + 1
 131 
 132     ABI_MALLOC(token,)
 133     token%id = n_symmetries
 134     call new_symmetry(token%data)
 135     token%next => my_symmetries
 136 
 137     my_symmetries => token
 138     if (AB_DBG) write(std_err,*) "AB symmetry: creation OK with id ", token%id
 139   end subroutine new_item
 140 
 141   subroutine free_item(token)
 142 
 143     type(symmetry_list), pointer :: token
 144 
 145     type(symmetry_list), pointer :: tmp
 146 
 147     if (.not. associated(token)) then
 148        return
 149     end if
 150 
 151     call free_symmetry(token%data)
 152 
 153     if (AB_DBG) write(std_err,*) "AB symmetry: free request on token ", token%id
 154     ! We remove token from the list.
 155     if (my_symmetries%id == token%id) then
 156        my_symmetries => token%next
 157     else
 158        tmp => my_symmetries
 159        do
 160           if (.not.associated(tmp)) then
 161              return
 162           end if
 163           if (associated(tmp%next) .and. tmp%next%id == token%id) then
 164              exit
 165           end if
 166           tmp => tmp%next
 167        end do
 168        tmp%next => token%next
 169     end if
 170     ABI_FREE(token)
 171     if (AB_DBG) write(std_err,*) "AB symmetry: free done"
 172   end subroutine free_item
 173 
 174   subroutine get_item(token, id)
 175 
 176     type(symmetry_list), pointer :: token
 177     integer, intent(in) :: id
 178 
 179     type(symmetry_list), pointer :: tmp
 180 
 181     if (AB_DBG) write(std_err,*) "AB symmetry: request list element ", id
 182     nullify(token)
 183 
 184     tmp => my_symmetries
 185     do
 186        if (.not. associated(tmp)) then
 187           exit
 188        end if
 189        if (tmp%id == id) then
 190           token => tmp
 191           return
 192        end if
 193        tmp => tmp%next
 194     end do
 195   end subroutine get_item
 196 
 197   subroutine symmetry_get_from_id(sym, id, errno)
 198 
 199     type(symmetry_type), pointer :: sym
 200     integer, intent(in) :: id
 201     integer, intent(out) :: errno
 202 
 203     type(symmetry_list), pointer :: token
 204 
 205     errno = AB7_NO_ERROR
 206     call get_item(token, id)
 207     if (associated(token)) then
 208        sym => token%data
 209        if (sym%nSym <= 0) then
 210           ! We do the computation of the matrix part.
 211           call compute_matrices(sym, errno)
 212        end if
 213     else
 214        errno = AB7_ERROR_OBJ
 215        nullify(sym)
 216     end if
 217   end subroutine symmetry_get_from_id
 218 
 219   subroutine new_symmetry(sym)
 220 
 221     type(symmetry_type), intent(out) :: sym
 222 
 223     if (AB_DBG) write(std_err,*) "AB symmetry: create a new symmetry object."
 224     nullify(sym%xRed)
 225     nullify(sym%spinAt)
 226     nullify(sym%typeAt)
 227     sym%tolsym   = tol8
 228     sym%auto     = .true.
 229     sym%nSym     = 0
 230     nullify(sym%sym)
 231     nullify(sym%symAfm)
 232     nullify(sym%transNon)
 233     sym%nBravSym = -1
 234     sym%withField   = .false.
 235     sym%withJellium = .false.
 236     sym%nzchempot = 0
 237     sym%withSpin = 1
 238     sym%withSpinOrbit = .false.
 239     sym%multiplicity = -1
 240     nullify(sym%indexingAtoms)
 241     sym%vacuum = 0
 242   end subroutine new_symmetry
 243 
 244   subroutine free_symmetry(sym)
 245 
 246     type(symmetry_type), intent(inout) :: sym
 247 
 248     if (AB_DBG) write(std_err,*) "AB symmetry: free a symmetry."
 249 
 250     if (associated(sym%xRed))  then
 251       ABI_FREE(sym%xRed)
 252     end if
 253     if (associated(sym%spinAt))  then
 254       ABI_FREE(sym%spinAt)
 255     end if
 256     if (associated(sym%typeAt))  then
 257       ABI_FREE(sym%typeAt)
 258     end if
 259     if (associated(sym%indexingAtoms))  then
 260       ABI_FREE(sym%indexingAtoms)
 261     end if
 262     if (associated(sym%sym))  then
 263       ABI_FREE(sym%sym)
 264     end if
 265     if (associated(sym%symAfm))  then
 266       ABI_FREE(sym%symAfm)
 267     end if
 268     if (associated(sym%transNon))  then
 269       ABI_FREE(sym%transNon)
 270     end if
 271   end subroutine free_symmetry
 272 
 273 
 274 
 275 
 276 
 277   subroutine symmetry_new(id)
 278 
 279     integer, intent(out) :: id
 280 
 281     type(symmetry_list), pointer :: token
 282 
 283     if (AB_DBG) write(std_err,*) "AB symmetry: call new symmetry."
 284     call new_item(token)
 285     id = token%id
 286   end subroutine symmetry_new
 287 
 288   subroutine symmetry_free(id)
 289 
 290     integer, intent(in) :: id
 291 
 292     type(symmetry_list), pointer :: token
 293 
 294     if (AB_DBG) write(std_err,*) "AB symmetry: call free symmetry."
 295 
 296     call get_item(token, id)
 297     if (associated(token)) then
 298       call free_item(token)
 299     end if
 300   end subroutine symmetry_free
 301 
 302   subroutine symmetry_set_tolerance(id, tolsym, errno)
 303 
 304     integer, intent(in) :: id
 305     real(dp), intent(in) :: tolsym
 306     integer, intent(out) :: errno
 307 
 308     type(symmetry_list), pointer :: token
 309 
 310     if (AB_DBG) write(std_err,*) "AB symmetry: call set tolerance."
 311 
 312     errno = AB7_NO_ERROR
 313     call get_item(token, id)
 314     if (.not. associated(token)) then
 315        errno = AB7_ERROR_OBJ
 316        return
 317     end if
 318 
 319     token%data%tolsym = tolsym
 320 
 321     ! We unset all the computed symmetries
 322     token%data%nBravSym = -1
 323     if (token%data%auto) then
 324        token%data%nSym  = 0
 325     end if
 326   end subroutine symmetry_set_tolerance
 327 
 328   subroutine symmetry_set_lattice(id, rprimd, errno)
 329 
 330     integer, intent(in) :: id
 331     real(dp), intent(in) :: rprimd(3,3)
 332     integer, intent(out) :: errno
 333 
 334     type(symmetry_list), pointer :: token
 335     real(dp) :: ucvol
 336     real(dp) :: gmet(3,3)
 337 
 338     if (AB_DBG) write(std_err,*) "AB symmetry: call set lattice."
 339     if (AB_DBG) write(std_err, "(A,3F12.6,A)") "  (", rprimd(:,1), ")"
 340     if (AB_DBG) write(std_err, "(A,3F12.6,A)") "  (", rprimd(:,2), ")"
 341     if (AB_DBG) write(std_err, "(A,3F12.6,A)") "  (", rprimd(:,3), ")"
 342 
 343     errno = AB7_NO_ERROR
 344     call get_item(token, id)
 345     if (.not. associated(token)) then
 346        errno = AB7_ERROR_OBJ
 347        return
 348     end if
 349 
 350     token%data%rprimd = rprimd
 351     call metric(gmet, token%data%gprimd, -1, token%data%rmet, rprimd, ucvol)
 352 
 353     ! We unset all the computed symmetries
 354     token%data%nBravSym = -1
 355     if (token%data%auto) then
 356        token%data%nSym  = 0
 357     end if
 358   end subroutine symmetry_set_lattice
 359 
 360   subroutine symmetry_set_structure(id, nAtoms, typeAt, xRed, errno)
 361 
 362     integer, intent(in) :: id
 363     integer, intent(in) :: nAtoms
 364     integer, intent(in) :: typeAt(nAtoms)
 365     real(dp), intent(in) :: xRed(3,nAtoms)
 366     integer, intent(out) :: errno
 367 
 368     type(symmetry_list), pointer :: token
 369     integer :: i
 370 
 371     if (AB_DBG) write(std_err,*) "AB symmetry: call set structure."
 372     if (AB_DBG) write(std_err, "(A,I3,A)") "  ", nAtoms, " atoms"
 373     if (AB_DBG) then
 374        do i = 1, nAtoms, 1
 375           if (AB_DBG) write(std_err,"(A,3F12.6,I3)") "  ", xRed(:, i), typeAt(i)
 376        end do
 377     end if
 378 
 379     errno = AB7_NO_ERROR
 380     call get_item(token, id)
 381     if (.not. associated(token)) then
 382        errno = AB7_ERROR_OBJ
 383        return
 384     end if
 385 
 386     token%data%nAtoms =  nAtoms
 387     ABI_MALLOC(token%data%typeAt,(nAtoms))
 388     token%data%typeAt = typeAt
 389     ABI_MALLOC(token%data%xRed,(3, nAtoms))
 390     token%data%xRed   = xRed
 391 
 392     ! We unset only the symmetries
 393     if (token%data%auto) then
 394        token%data%nSym = 0
 395     end if
 396     if (associated(token%data%indexingAtoms)) then
 397       ABI_FREE(token%data%indexingAtoms)
 398     end if
 399   end subroutine symmetry_set_structure
 400 
 401   subroutine symmetry_set_spin(id, nAtoms, spinAt, errno)
 402 
 403     integer, intent(in) :: id
 404     integer, intent(in) :: nAtoms
 405     real(dp), intent(in) :: spinAt(3,nAtoms)
 406     integer, intent(out) :: errno
 407 
 408     type(symmetry_list), pointer :: token
 409     integer :: i
 410 
 411     if (AB_DBG) write(std_err,*) "AB symmetry: call set spin."
 412     if (AB_DBG) then
 413        do i = 1, nAtoms, 1
 414           if (AB_DBG) write(std_err,"(A,3F12.6)") "  ", spinAt(:, i)
 415        end do
 416     end if
 417 
 418     errno = AB7_NO_ERROR
 419     call get_item(token, id)
 420     if (.not. associated(token)) then
 421        errno = AB7_ERROR_OBJ
 422        return
 423     end if
 424     if (token%data%nAtoms /= nAtoms) then
 425        errno = AB7_ERROR_ARG
 426        return
 427     end if
 428 
 429     token%data%withSpin = 4
 430     ABI_MALLOC(token%data%spinAt,(3, nAtoms))
 431     token%data%spinAt = spinAt
 432 
 433     ! We unset only the symmetries
 434     if (token%data%auto) then
 435        token%data%nSym  = 0
 436     end if
 437   end subroutine symmetry_set_spin
 438 
 439   subroutine symmetry_set_collinear_spin(id, nAtoms, spinAt, errno)
 440 
 441     integer, intent(in) :: id
 442     integer, intent(in) :: nAtoms
 443     integer, intent(in) :: spinAt(nAtoms)
 444     integer, intent(out) :: errno
 445 
 446     type(symmetry_list), pointer :: token
 447     integer :: i
 448 
 449     if (AB_DBG) write(std_err,*) "AB symmetry: call set collinear spin."
 450     if (AB_DBG) then
 451        do i = 1, nAtoms, 1
 452           if (AB_DBG) write(std_err,"(A,I3)") "  ", spinAt(i)
 453        end do
 454     end if
 455 
 456     errno = AB7_NO_ERROR
 457     call get_item(token, id)
 458     if (.not. associated(token)) then
 459        errno = AB7_ERROR_OBJ
 460        return
 461     end if
 462     if (token%data%nAtoms /= nAtoms) then
 463        errno = AB7_ERROR_ARG
 464        return
 465     end if
 466 
 467     token%data%withSpin = 2
 468     ABI_MALLOC(token%data%spinAt,(3,nAtoms))
 469     token%data%spinAt = 0._dp
 470     token%data%spinAt(3, :) = real(spinAt)
 471 
 472     ! We unset only the symmetries
 473     if (token%data%auto) then
 474        token%data%nSym  = 0
 475     end if
 476   end subroutine symmetry_set_collinear_spin
 477 
 478   subroutine symmetry_set_spin_orbit(id, withSpinOrbit, errno)
 479 
 480     integer, intent(in) :: id
 481     logical, intent(in) :: withSpinOrbit
 482     integer, intent(out) :: errno
 483 
 484     type(symmetry_list), pointer :: token
 485 
 486     if (AB_DBG) write(std_err,*) "AB symmetry: call set spin orbit."
 487 
 488     errno = AB7_NO_ERROR
 489     call get_item(token, id)
 490     if (.not. associated(token)) then
 491        errno = AB7_ERROR_OBJ
 492        return
 493     end if
 494 
 495     token%data%withSpinOrbit = withSpinOrbit
 496 
 497     ! We unset only the symmetries
 498     if (token%data%auto) then
 499        token%data%nSym  = 0
 500     end if
 501   end subroutine symmetry_set_spin_orbit
 502 
 503   subroutine symmetry_set_field(id, field, errno)
 504 
 505     integer, intent(in) :: id
 506     real(dp), intent(in) :: field(3)
 507     integer, intent(out) :: errno
 508 
 509     type(symmetry_list), pointer :: token
 510 
 511     if (AB_DBG) write(std_err,*) "AB symmetry: call set field."
 512 
 513     errno = AB7_NO_ERROR
 514     call get_item(token, id)
 515     if (.not. associated(token)) then
 516        errno = AB7_ERROR_OBJ
 517        return
 518     end if
 519 
 520     token%data%withField = .true.
 521     token%data%field = field
 522 
 523     ! We unset all the computed symmetries
 524     token%data%nBravSym = -1
 525     if (token%data%auto) then
 526        token%data%nSym  = 0
 527     end if
 528   end subroutine symmetry_set_field
 529 
 530   subroutine symmetry_set_jellium(id, jellium, errno)
 531 
 532     integer, intent(in) :: id
 533     logical, intent(in) :: jellium
 534     integer, intent(out) :: errno
 535 
 536     type(symmetry_list), pointer :: token
 537 
 538     if (AB_DBG) write(std_err,*) "AB symmetry: call set jellium."
 539 
 540     errno = AB7_NO_ERROR
 541     call get_item(token, id)
 542     if (.not. associated(token)) then
 543        errno = AB7_ERROR_OBJ
 544        return
 545     end if
 546 
 547     token%data%withJellium = jellium
 548 
 549     ! We unset only the symmetries
 550     if (token%data%auto) then
 551        token%data%nSym  = 0
 552     end if
 553   end subroutine symmetry_set_jellium
 554 
 555 
 556   subroutine symmetry_set_nzchempot(id, nzchempot, errno)
 557 
 558     integer, intent(in) :: id
 559     integer, intent(in) :: nzchempot
 560     integer, intent(out) :: errno
 561 
 562     type(symmetry_list), pointer :: token
 563 
 564     if (AB_DBG) write(std_err,*) "AB symmetry: call set nzchempot."
 565 
 566     errno = AB7_NO_ERROR
 567     call get_item(token, id)
 568     if (.not. associated(token)) then
 569        errno = AB7_ERROR_OBJ
 570        return
 571     end if
 572 
 573     token%data%nzchempot = nzchempot
 574 
 575     ! We unset only the symmetries
 576     if (token%data%auto) then
 577        token%data%nSym  = 0
 578     end if
 579   end subroutine symmetry_set_nzchempot
 580 
 581   subroutine symmetry_set_periodicity(id, periodic, errno)
 582 
 583     integer, intent(in) :: id
 584     logical, intent(in) :: periodic(3)
 585     integer, intent(out) :: errno
 586 
 587     type(symmetry_list), pointer :: token
 588 
 589     if (AB_DBG) write(std_err,*) "AB symmetry: call set periodicity."
 590     if (AB_DBG) write(std_err, "(A,3L1,A)") "  (", periodic, ")"
 591 
 592     errno = AB7_NO_ERROR
 593     call get_item(token, id)
 594     if (.not. associated(token)) then
 595        errno = AB7_ERROR_OBJ
 596        return
 597     end if
 598 
 599     token%data%vacuum = 0
 600     if (.not. periodic(1)) token%data%vacuum(1) = 1
 601     if (.not. periodic(2)) token%data%vacuum(2) = 1
 602     if (.not. periodic(3)) token%data%vacuum(3) = 1
 603   end subroutine symmetry_set_periodicity
 604 
 605 
 606 
 607 
 608 
 609   subroutine symmetry_get_n_atoms(id, nAtoms, errno)
 610     !scalars
 611     integer, intent(in) :: id
 612     integer, intent(out) :: errno
 613     integer, intent(out) :: nAtoms
 614 
 615     type(symmetry_list), pointer :: token
 616 
 617     if (AB_DBG) write(std_err,*) "AB symmetry: call get nAtoms."
 618 
 619     errno = AB7_NO_ERROR
 620     call get_item(token, id)
 621     if (.not. associated(token)) then
 622        errno = AB7_ERROR_OBJ
 623        return
 624     end if
 625 
 626     nAtoms = token%data%nAtoms
 627   end subroutine symmetry_get_n_atoms
 628 
 629   subroutine compute_bravais(sym)
 630 
 631     type(symmetry_type), intent(inout) :: sym
 632 
 633     integer :: berryopt
 634 
 635     ! We do the computation
 636     if (sym%withField) then
 637        berryopt = 4
 638     else
 639        berryopt = 0
 640     end if
 641     if (AB_DBG) write(std_err,*) "AB symmetry: call ABINIT symlatt."
 642     call symlatt(sym%bravais, std_out, AB7_MAX_SYMMETRIES, &
 643          & sym%nBravSym, sym%bravSym, sym%rprimd, sym%tolsym)
 644     if (AB_DBG) write(std_err,*) "AB symmetry: call ABINIT OK."
 645     if (AB_DBG) write(std_err, "(A,I3)") "  nSymBrav :", sym%nBravSym
 646     if (AB_DBG) write(std_err, "(A,I3)") "  holohedry:", sym%bravais(1)
 647     if (AB_DBG) write(std_err, "(A,I3)") "  center   :", sym%bravais(2)
 648   end subroutine compute_bravais
 649 
 650   subroutine symmetry_get_bravais(id, bravais, holohedry, center, &
 651        & nBravSym, bravSym, errno)
 652     !scalars
 653     integer, intent(in) :: id
 654     integer, intent(out) :: errno
 655     integer, intent(out) :: nBravSym, holohedry, center
 656     !arrays
 657     integer, intent(out) :: bravais(3,3), bravSym(3, 3, AB7_MAX_SYMMETRIES)
 658 
 659     type(symmetry_list), pointer :: token
 660 
 661     if (AB_DBG) write(std_err,*) "AB symmetry: call get bravais."
 662 
 663     errno = AB7_NO_ERROR
 664     call get_item(token, id)
 665     if (.not. associated(token)) then
 666        errno = AB7_ERROR_OBJ
 667        return
 668     end if
 669 
 670     if (token%data%nBravSym < 0) then
 671        ! We do the computation
 672        call compute_bravais(token%data)
 673     end if
 674 
 675     holohedry = token%data%bravais(1)
 676     center    = token%data%bravais(2)
 677     bravais   = reshape(token%data%bravais(3:11), (/ 3,3 /))
 678     nBravSym  = token%data%nBravSym
 679     bravSym(:, :, 1:nBravSym) = token%data%bravSym(:, :, 1:nBravSym)
 680   end subroutine symmetry_get_bravais
 681 
 682   subroutine compute_matrices(sym, errno)
 683 
 684     type(symmetry_type), intent(inout) :: sym
 685     integer, intent(out) :: errno
 686 
 687     integer :: berryopt, invar_z
 688     integer :: use_inversion
 689     real(dp), pointer :: spinAt_(:,:)
 690     integer  :: sym_(3, 3, AB7_MAX_SYMMETRIES)
 691     real(dp) :: transNon_(3, AB7_MAX_SYMMETRIES)
 692     integer  :: symAfm_(AB7_MAX_SYMMETRIES)
 693 
 694     errno = AB7_NO_ERROR
 695 
 696     if (sym%nBravSym < 0) then
 697        ! We do the computation of the Bravais part.
 698        call compute_bravais(sym)
 699     end if
 700 
 701     if (sym%withField) then
 702        berryopt = 4
 703     else
 704        berryopt = 0
 705     end if
 706     if (sym%withJellium .or. sym%nzchempot/=0) then
 707        invar_z = 2
 708     else
 709        invar_z = 0
 710     end if
 711     if (sym%withSpin == 2 .or. sym%withSpin == 4) then
 712        spinAt_ => sym%spinAt
 713     else
 714        ABI_MALLOC(spinAt_,(3, sym%nAtoms))
 715        spinAt_ = 0
 716     end if
 717     if (sym%withSpinOrbit) then
 718        use_inversion = 0
 719     else
 720        use_inversion = 1
 721     end if
 722 
 723     if (sym%nsym == 0) then
 724        if (AB_DBG) write(std_err,*) "AB symmetry: call ABINIT symfind."
 725        call symfind(sym%gprimd, AB7_MAX_SYMMETRIES, &
 726             & sym%nAtoms, sym%nBravSym, sym%withSpin, sym%nSym, 0, sym%bravSym, spinAt_, &
 727             & symAfm_, sym_, transNon_, sym%tolsym, sym%typeAt, &
 728             & use_inversion, sym%xRed, invardir_red=sym%field, invar_z=invar_z)
 729        if (AB_DBG) write(std_err,*) "AB symmetry: call ABINIT OK."
 730        if (AB_DBG) write(std_err, "(A,I3)") "  nSym:", sym%nSym
 731        if (associated(sym%sym))  then
 732           ABI_FREE(sym%sym)
 733        end if
 734        if (associated(sym%symAfm))  then
 735           ABI_FREE(sym%symAfm)
 736        end if
 737        if (associated(sym%transNon))  then
 738           ABI_FREE(sym%transNon)
 739        end if
 740        ABI_MALLOC(sym%sym,(3, 3, sym%nSym))
 741        sym%sym(:,:,:) = sym_(:,:, 1:sym%nSym)
 742        ABI_MALLOC(sym%symAfm,(sym%nSym))
 743        sym%symAfm(:) = symAfm_(1:sym%nSym)
 744        ABI_MALLOC(sym%transNon,(3, sym%nSym))
 745        sym%transNon(:,:) = transNon_(:, 1:sym%nSym)
 746     else if (sym%nsym < 0) then
 747        sym%nsym = -sym%nsym
 748        sym_(:,:, 1:sym%nSym) = sym%sym(:,:,:)
 749        transNon_(:, 1:sym%nSym) = sym%transNon(:,:)
 750        symAfm_(1:sym%nSym) = sym%symAfm(:)
 751     end if
 752 
 753     if (sym%withSpin == 1) then
 754        ABI_FREE(spinAt_)
 755     end if
 756 
 757     if (AB_DBG) write(std_err,*) "AB symmetry: call ABINIT symanal."
 758     call symanal(sym%bravais, 0, sym%genAfm, AB7_MAX_SYMMETRIES, sym%nSym, &
 759          & sym%pointGroupMagn, sym%rprimd, sym%spaceGroup, symAfm_, &
 760          & sym_, transNon_, sym%tolsym)
 761     if (AB_DBG) write(std_err,*) "AB symmetry: call ABINIT OK."
 762     sym%transNon(:,:) = transNon_(:, 1:sym%nSym)
 763 
 764     if (sym%bravais(1) < 0) then
 765        sym%multiplicity = 2
 766     else
 767        sym%multiplicity = 1
 768     end if
 769     if (AB_DBG) write(std_err, "(A,I3)") "  multi:", sym%multiplicity
 770     if (AB_DBG) write(std_err, "(A,I3)") "  space:", sym%spaceGroup
 771   end subroutine compute_matrices
 772 
 773   subroutine symmetry_get_n_sym(id, nSym, errno)
 774     !scalars
 775     integer, intent(in) :: id
 776     integer, intent(out) :: errno
 777     integer, intent(out) :: nSym
 778 
 779     type(symmetry_list), pointer :: token
 780 
 781     if (AB_DBG) write(std_err,*) "AB symmetry: call get nSym."
 782 
 783     errno = AB7_NO_ERROR
 784     call get_item(token, id)
 785     if (.not. associated(token)) then
 786        errno = AB7_ERROR_OBJ
 787        return
 788     end if
 789 
 790     if (token%data%nSym <= 0) then
 791        ! We do the computation of the matrix part.
 792        call compute_matrices(token%data, errno)
 793     end if
 794 
 795     nSym = token%data%nSym
 796   end subroutine symmetry_get_n_sym
 797 
 798   subroutine symmetry_set_n_sym(id, nSym, sym, transNon, symAfm, errno)
 799     !scalars
 800     integer, intent(in)  :: id
 801     integer, intent(in)  :: nSym
 802     integer, intent(in)  :: sym(3, 3, nSym)
 803     real(dp), intent(in) :: transNon(3, nSym)
 804     integer, intent(in)  :: symAfm(nSym)
 805     integer, intent(out) :: errno
 806 
 807     type(symmetry_list), pointer :: token
 808 
 809     if (AB_DBG) write(std_err,*) "AB symmetry: call get nSym."
 810 
 811     errno = AB7_NO_ERROR
 812     call get_item(token, id)
 813     if (.not. associated(token)) then
 814        errno = AB7_ERROR_OBJ
 815        return
 816     end if
 817 
 818     if (nSym <= 0) then
 819        errno = AB7_ERROR_ARG
 820        return
 821     else
 822        ABI_MALLOC(token%data%sym, (3, 3, nSym))
 823        token%data%sym(:,:,:) = sym(:,:,:)
 824        ABI_MALLOC(token%data%symAfm, (nSym))
 825        token%data%symAfm(:) = symAfm(:)
 826        ABI_MALLOC(token%data%transNon, (3, nSym))
 827        token%data%transNon(:,:) = transNon(:,:)
 828 
 829        token%data%auto = .false.
 830        token%data%nsym = -nSym
 831     end if
 832 
 833     ! We do the computation of the matrix part.
 834     call compute_matrices(token%data, errno)
 835   end subroutine symmetry_set_n_sym
 836 
 837   subroutine symmetry_get_matrices(id, nSym, sym, transNon, symAfm, errno)
 838 
 839     integer, intent(in) :: id
 840     integer, intent(out) :: errno
 841     integer, intent(out) :: nSym
 842     integer, intent(out)  :: sym(3, 3, AB7_MAX_SYMMETRIES)
 843     integer, intent(out)  :: symAfm(AB7_MAX_SYMMETRIES)
 844     real(dp), intent(out) :: transNon(3, AB7_MAX_SYMMETRIES)
 845 
 846     type(symmetry_list), pointer :: token
 847 
 848     if (AB_DBG) write(std_err,*) "AB symmetry: call get matrices."
 849 
 850     errno = AB7_NO_ERROR
 851     call get_item(token, id)
 852     if (.not. associated(token)) then
 853        errno = AB7_ERROR_OBJ
 854        return
 855     end if
 856 
 857     if (token%data%nSym <= 0) then
 858        ! We do the computation of the matrix part.
 859        call compute_matrices(token%data, errno)
 860     end if
 861 
 862     nSym                = token%data%nSym
 863     sym(:, :, 1:nSym)   = token%data%sym(:, :,:)
 864     symAfm(1:nSym)      = token%data%symAfm(:)
 865     transNon(:, 1:nSym) = token%data%transNon(:,:)
 866   end subroutine symmetry_get_matrices
 867 
 868   subroutine symmetry_get_matrices_p(id, nSym, sym, transNon, symAfm, errno)
 869 
 870     integer, intent(in) :: id
 871     integer, intent(out) :: errno
 872     integer, intent(out) :: nSym
 873     integer, pointer  :: sym(:,:,:)
 874     integer, pointer  :: symAfm(:)
 875     real(dp), pointer :: transNon(:,:)
 876 
 877     type(symmetry_list), pointer :: token
 878 
 879     if (AB_DBG) write(std_err,*) "AB symmetry: call get matrices as pointers."
 880 
 881     errno = AB7_NO_ERROR
 882     call get_item(token, id)
 883     if (.not. associated(token)) then
 884        errno = AB7_ERROR_OBJ
 885        return
 886     end if
 887 
 888     if (token%data%nSym <= 0) then
 889        ! We do the computation of the matrix part.
 890        call compute_matrices(token%data, errno)
 891     end if
 892 
 893     nSym     =  token%data%nSym
 894     sym      => token%data%sym
 895     symAfm   => token%data%symAfm
 896     transNon => token%data%transNon
 897   end subroutine symmetry_get_matrices_p
 898 
 899   subroutine symmetry_get_multiplicity(id, multiplicity, errno)
 900 
 901     integer, intent(in) :: id
 902     integer, intent(out) :: multiplicity, errno
 903 
 904     type(symmetry_list), pointer :: token
 905 
 906     if (AB_DBG) write(std_err,*) "AB symmetry: call get multiplicity."
 907 
 908     errno = AB7_NO_ERROR
 909     call get_item(token, id)
 910     if (.not. associated(token)) then
 911        errno = AB7_ERROR_OBJ
 912        return
 913     end if
 914 
 915     if (token%data%multiplicity < 0) then
 916        ! We do the computation of the matrix part.
 917        call compute_matrices(token%data, errno)
 918     end if
 919     multiplicity = token%data%multiplicity
 920   end subroutine symmetry_get_multiplicity
 921 
 922   subroutine symmetry_get_group(id, spaceGroup, spaceGroupId, &
 923        & pointGroupMagn, genAfm, errno)
 924 
 925     integer, intent(in)            :: id
 926     integer, intent(out)           :: errno
 927     real(dp), intent(out)          :: genAfm(3)
 928     character(len=15), intent(out) :: spaceGroup
 929     integer, intent(out)           :: spaceGroupId, pointGroupMagn
 930 
 931     type(symmetry_list), pointer  :: token
 932     integer :: sporder
 933     character(len=1)  :: brvLattice
 934     character(len=15) :: ptintsb,ptschsb,schsb,spgrp
 935     character(len=35) :: intsbl
 936 
 937     if (AB_DBG) write(std_err,*) "AB symmetry: call get group."
 938 
 939     errno = AB7_NO_ERROR
 940     call get_item(token, id)
 941     if (.not. associated(token)) then
 942        errno = AB7_ERROR_OBJ
 943        return
 944     end if
 945 
 946     if (token%data%multiplicity < 0) then
 947        ! We do the computation of the matrix part.
 948        call compute_matrices(token%data, errno)
 949     end if
 950 
 951     if (token%data%multiplicity /= 1) then
 952        errno = AB7_ERROR_SYM_NOT_PRIMITIVE
 953        return
 954     end if
 955 
 956     call spgdata(brvLattice,spgrp,intsbl,ptintsb,ptschsb,&
 957          &  schsb,1,token%data%spaceGroup,sporder,1)
 958 
 959     write(spaceGroup, "(3A)") brvLattice, " ", trim(spgrp(1:13))
 960     pointGroupMagn = token%data%pointGroupMagn
 961     spaceGroupId   = token%data%spaceGroup
 962     genAfm         = token%data%genAfm
 963   end subroutine symmetry_get_group
 964 
 965   subroutine compute_equivalent_atoms(sym)
 966 
 967     type(symmetry_type), intent(inout) :: sym
 968 
 969     integer, allocatable :: symrec(:,:,:)
 970     integer :: isym
 971 
 972     if (.not. associated(sym%indexingAtoms)) then
 973       ABI_MALLOC(sym%indexingAtoms,(4, sym%nSym, sym%nAtoms))
 974     end if
 975 
 976     !Get the symmetry matrices in terms of reciprocal basis
 977     ABI_MALLOC(symrec,(3, 3, sym%nSym))
 978     do isym = 1, sym%nSym, 1
 979        call mati3inv(sym%sym(:,:,isym), symrec(:,:,isym))
 980     end do
 981 
 982     !Obtain a list of rotated atom labels:
 983     call symatm(sym%indexingAtoms, sym%nAtoms, sym%nSym, symrec, &
 984          & sym%transNon, sym%tolsym, sym%typeAt, sym%xRed)
 985 
 986     ABI_FREE(symrec)
 987   end subroutine compute_equivalent_atoms
 988 
 989   subroutine symmetry_get_equivalent_atom(id, equiv, iAtom, errno)
 990 
 991     integer, intent(in)  :: id
 992     integer, intent(in)  :: iAtom
 993     integer, intent(out) :: equiv(4, AB7_MAX_SYMMETRIES)
 994     integer, intent(out) :: errno
 995 
 996     type(symmetry_list), pointer  :: token
 997 
 998     if (AB_DBG) write(std_err,*) "AB symmetry: call get equivalent."
 999 
1000     errno = AB7_NO_ERROR
1001     call get_item(token, id)
1002     if (.not. associated(token)) then
1003        errno = AB7_ERROR_OBJ
1004        return
1005     end if
1006 
1007     if (iAtom < 1 .or. iAtom > token%data%nAtoms) then
1008        errno = AB7_ERROR_ARG
1009        return
1010     end if
1011 
1012     if (.not. associated(token%data%indexingAtoms)) then
1013        ! We do the computation of the matrix part.
1014        call compute_equivalent_atoms(token%data)
1015     end if
1016 
1017     equiv(:, 1:token%data%nSym) = token%data%indexingAtoms(:,:,iAtom)
1018   end subroutine symmetry_get_equivalent_atom
1019 
1020   subroutine symmetry_get_type(id, iSym, label, type, errno)
1021 
1022     integer, intent(in)  :: id, iSym
1023     character(len = 128), intent(out) :: label
1024     integer, intent(out) :: errno, type
1025 
1026     integer :: det
1027     type(symmetry_list), pointer :: token
1028     type(symmetry_type), pointer :: sym
1029 
1030     if (AB_DBG) write(std_err,*) "AB symmetry: call get type."
1031 
1032     errno = AB7_NO_ERROR
1033     call get_item(token, id)
1034     if (.not. associated(token)) then
1035        errno = AB7_ERROR_OBJ
1036        return
1037     end if
1038     sym => token%data
1039 
1040     if (iSym < 1 .or. iSym > sym%nSym) then
1041        errno = AB7_ERROR_ARG
1042        return
1043     end if
1044 
1045     if (sym%multiplicity < 0) then
1046        ! We do the computation of the matrix part.
1047        call compute_matrices(sym, errno)
1048     end if
1049 
1050     ! Calculate determinant.
1051     call mati3det(sym%sym(:,:,iSym),det)
1052     call symcharac(sym%bravais(2), det, sym%bravais(1), iSym, label, &
1053          sym%sym(:,:,iSym), sym%transNon(:, iSym), type)
1054   end subroutine symmetry_get_type
1055 
1056 end module m_ab7_symmetry