The Elk Code
bectask.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine bectask(fnum,fext)
7 use modmain
8 use modphonon
9 use modmpi
10 implicit none
11 ! arguments
12 integer, intent(in) :: fnum
13 character(*), intent(out) :: fext
14 ! local variables
15 logical exist
16 ! only master process should search for file
17 if (.not.mp_mpi) goto 10
18 do ipph=1,3
19  do isph=1,nspecies
20  do iaph=1,natoms(isph)
21 ! Born effective charge file extension
22  call becfext(isph,iaph,ipph,fext)
23 ! determine if the BEC file with this extension exists
24  inquire(file='BEC'//trim(fext),exist=exist)
25  if (.not.exist) then
26  open(fnum,file='BEC'//trim(fext),form='FORMATTED')
28  goto 10
29  end if
30  end do
31  end do
32 end do
33 isph=0; iaph=0; iasph=0; ipph=0
34 write(*,'("Info(bectask): nothing more to do")')
35 10 continue
36 ! broadcast to all other MPI processes
37 call mpi_bcast(isph,1,mpi_integer,0,mpicom,ierror)
38 call mpi_bcast(iaph,1,mpi_integer,0,mpicom,ierror)
39 call mpi_bcast(iasph,1,mpi_integer,0,mpicom,ierror)
40 call mpi_bcast(ipph,1,mpi_integer,0,mpicom,ierror)
41 if (isph == 0) then
42  fext='.OUT'
43 else
44  call becfext(isph,iaph,ipph,fext)
45 end if
46 end subroutine
47 
logical mp_mpi
Definition: modmpi.f90:17
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
integer iasph
Definition: modphonon.f90:15
integer ipph
Definition: modphonon.f90:15
integer isph
Definition: modphonon.f90:15
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
Definition: modmpi.f90:6
integer nspecies
Definition: modmain.f90:34
subroutine bectask(fnum, fext)
Definition: bectask.f90:7
integer iaph
Definition: modphonon.f90:15
subroutine becfext(is, ia, ip, fext)
Definition: becfext.f90:7
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19