TABLE OF CONTENTS


defs_wvltypes/wvl_descr_atoms_set [ Functions ]

[ Top ] [ defs_wvltypes ] [ Functions ]

NAME

 wvl_descr_atoms_set

FUNCTION

 Defines wvl%atoms% data structure

COPYRIGHT

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

INPUTS

 acell(3)=unit cell length scales (bohr)
 dtset <type(dataset_type)>=all input variables for this dataset

OUTPUT

 wvl <type(wvl_internal_type)>= wavelet type
                 | nat      =  number of atoms
                 | ntypes   =  number of species
                 | alat1    =  acell(1)
                 | alat2    =  acell(2)
                 | alat3    =  acell(3)
                 | iatype   =  types for atoms
                 | lfrztyp  =  flag for the movement of atoms.
                 | natpol   =  integer related to polarisation at the first step

PARENTS

      gstate,wvl_memory

CHILDREN

      allocate_atoms_nat,allocate_atoms_ntypes,astruct_set_n_atoms
      astruct_set_n_types,f_release_routine,f_routine

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 subroutine wvl_descr_atoms_set(acell, icoulomb, natom, ntypat, typat, wvl)
 46 
 47  use m_profiling_abi
 48  use m_errors
 49 
 50  use defs_basis
 51  use defs_wvltypes
 52 #if defined HAVE_BIGDFT
 53  use BigDFT_API, only: atoms_data_null,f_routine,f_release_routine,&
 54 &                      astruct_set_n_atoms,astruct_set_n_types,&
 55 &                      allocate_atoms_nat,allocate_atoms_ntypes
 56 #endif
 57 
 58 !This section has been created automatically by the script Abilint (TD).
 59 !Do not modify the following lines by hand.
 60 #undef ABI_FUNC
 61 #define ABI_FUNC 'wvl_descr_atoms_set'
 62 !End of the abilint section
 63 
 64   implicit none
 65 
 66 !Arguments ------------------------------------
 67 !scalars
 68   integer, intent(in)                    :: icoulomb, natom, ntypat
 69   type(wvl_internal_type), intent(inout) :: wvl
 70   !arrays
 71   integer, intent(in)                    :: typat(natom)
 72   real(dp), intent(in)                   :: acell(3)
 73 
 74 !Local variables-------------------------------
 75 !scalars
 76 #if defined HAVE_BIGDFT
 77   integer :: itype
 78 #endif
 79 
 80 ! *********************************************************************
 81 
 82 #if defined HAVE_BIGDFT
 83 
 84  call f_routine(ABI_FUNC)
 85 
 86  wvl%atoms=atoms_data_null()
 87 
 88 !We create the atoms_data structure from this dataset
 89 !to be used later in BigDFT routines.
 90  if (icoulomb == 0) then
 91    wvl%atoms%astruct%geocode = 'P'
 92  else if (icoulomb == 1) then
 93    wvl%atoms%astruct%geocode = 'F'
 94  else if (icoulomb == 2) then
 95    wvl%atoms%astruct%geocode = 'S'
 96  end if
 97  write(wvl%atoms%astruct%units, "(A)") "Bohr"
 98 
 99  call astruct_set_n_atoms(wvl%atoms%astruct, natom)
100  call astruct_set_n_types(wvl%atoms%astruct, ntypat)
101 
102  do itype = 1, ntypat, 1
103    write(wvl%atoms%astruct%atomnames(itype), "(A,I2)") "At. type", itype
104  end do
105  wvl%atoms%astruct%cell_dim(1)   =  acell(1)
106  wvl%atoms%astruct%cell_dim(2)   =  acell(2)
107  wvl%atoms%astruct%cell_dim(3)   =  acell(3)
108  wvl%atoms%astruct%iatype   = typat
109 
110  wvl%atoms%astruct%sym%symObj = 0
111 
112  call allocate_atoms_nat(wvl%atoms)
113  call allocate_atoms_ntypes(wvl%atoms)
114 
115  call f_release_routine()
116 
117 #else
118  BIGDFT_NOTENABLED_ERROR()
119  if (.false.) write(std_out,*) icoulomb,natom,ntypat,wvl%h(1),typat(1),acell(1)
120 #endif
121  
122 end subroutine wvl_descr_atoms_set