File: dyntask.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (58 lines) | stat: -rw-r--r-- 1,524 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine dyntask(fnum,fext)
use modmain
use modphonon
use modmpi
implicit none
! arguments
integer, intent(in) :: fnum
character(*), intent(out) :: fext
! local variables
logical exist
! only master process should search for file
if (.not.mp_mpi) goto 10
do ipph=1,3
  do isph=1,nspecies
    do iaph=1,natoms(isph)
      do iqph=1,nqpt
! construct the phonon file extension
        call phfext(iqph,isph,iaph,ipph,fext)
! determine if the DYN file with this extension exists
        inquire(file='DYN'//trim(fext),exist=exist)
        if (.not.exist) then
          open(fnum,file='DYN'//trim(fext),form='FORMATTED')
          iasph=idxas(iaph,isph)
          goto 10
        end if
      end do
    end do
  end do
end do
iqph=0; isph=0; iaph=0; iasph=0; ipph=0
write(*,*)
write(*,'("Info(dyntask): nothing more to do")')
10 continue
! broadcast to all other MPI processes
call mpi_bcast(iqph,1,mpi_integer,0,mpicom,ierror)
call mpi_bcast(isph,1,mpi_integer,0,mpicom,ierror)
call mpi_bcast(iaph,1,mpi_integer,0,mpicom,ierror)
call mpi_bcast(iasph,1,mpi_integer,0,mpicom,ierror)
call mpi_bcast(ipph,1,mpi_integer,0,mpicom,ierror)
if (iqph.eq.0) then
  fext='.OUT'
else
  call phfext(iqph,isph,iaph,ipph,fext)
end if
! set the q=0 flag
if (iqph.eq.iq0) then
  tphq0=.true.
else
  tphq0=.false.
end if
return
end subroutine