TABLE OF CONTENTS


ABINIT/pspheads_comm [ Functions ]

[ Top ] [ Functions ]

NAME

 pspheads_comm

FUNCTION

 Communicate pspheads to all processors

COPYRIGHT

 Copyright (C) 2009-2018 ABINIT group (DCA, XG, GMR, FrD, AF, MT)
 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

  npsp=number of pseudopotentials
  test_paw=0 if no PAW, 1 if PAW

SIDE EFFECTS

  pspheads(npsp)=<type pspheader_type>=all the important information from the
   pseudopotential file headers, as well as the psp file names. On one processor at input,
   on all processors at output

PARENTS

      m_ab7_invars_f90

CHILDREN

      timab,xmpi_bcast

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 subroutine pspheads_comm(npsp,pspheads,test_paw)
 40 
 41  use defs_basis
 42  use defs_datatypes
 43  use m_profiling_abi
 44  use m_xmpi
 45 
 46 !This section has been created automatically by the script Abilint (TD).
 47 !Do not modify the following lines by hand.
 48 #undef ABI_FUNC
 49 #define ABI_FUNC 'pspheads_comm'
 50  use interfaces_18_timing
 51 !End of the abilint section
 52 
 53  implicit none
 54 
 55 !Arguments ------------------------------------
 56  integer,intent(in) :: npsp
 57  integer,intent(inout) :: test_paw
 58  type(pspheader_type),intent(inout) :: pspheads(npsp)
 59 
 60 !Local variables-------------------------------
 61 #if defined HAVE_MPI
 62 !scalars
 63  integer,parameter :: master=0
 64  integer :: ierr,comm
 65 !arrays
 66  integer,allocatable :: list_int(:)
 67  real(dp) :: tsec(2)
 68  real(dp),allocatable :: list_dpr(:)
 69  character(len=fnlen),allocatable :: list_char(:)
 70 #endif
 71 
 72 !*************************************************************************
 73 
 74 #if defined HAVE_MPI
 75  call timab(48,1,tsec)
 76 
 77  comm = xmpi_world
 78 
 79 !Broadcast the characters (file names and titles)
 80  ABI_ALLOCATE(list_char,(3*npsp))
 81  list_char(1:npsp)=pspheads(1:npsp)%filpsp
 82  list_char(npsp+1:2*npsp)=pspheads(1:npsp)%title
 83  list_char(2*npsp+1:3*npsp)=pspheads(1:npsp)%md5_checksum
 84 
 85  call xmpi_bcast(list_char,master,comm,ierr)
 86 
 87  pspheads(1:npsp)%filpsp=list_char(1:npsp)
 88  pspheads(1:npsp)%title=list_char(npsp+1:2*npsp)
 89  pspheads(1:npsp)%md5_checksum=list_char(2*npsp+1:3*npsp)(1:md5_slen)
 90  ABI_DEALLOCATE(list_char)
 91 
 92 !Brodcast the integers
 93  ABI_ALLOCATE(list_int,(1+13*npsp))
 94  list_int(1        :   npsp) = pspheads(1:npsp)%nproj(0)
 95  list_int(1+   npsp: 2*npsp) = pspheads(1:npsp)%nproj(1)
 96  list_int(1+ 2*npsp: 3*npsp) = pspheads(1:npsp)%nproj(2)
 97  list_int(1+ 3*npsp: 4*npsp) = pspheads(1:npsp)%nproj(3)
 98  list_int(1+ 4*npsp: 5*npsp) = pspheads(1:npsp)%lmax
 99  list_int(1+ 5*npsp: 6*npsp) = pspheads(1:npsp)%xccc
100  list_int(1+ 6*npsp: 7*npsp) = pspheads(1:npsp)%pspxc
101  list_int(1+ 7*npsp: 8*npsp) = pspheads(1:npsp)%pspdat
102  list_int(1+ 8*npsp: 9*npsp) = pspheads(1:npsp)%pspcod
103  list_int(1+ 9*npsp:10*npsp) = pspheads(1:npsp)%pspso
104  list_int(1+10*npsp:11*npsp) = pspheads(1:npsp)%nprojso(1)
105  list_int(1+11*npsp:12*npsp) = pspheads(1:npsp)%nprojso(2)
106  list_int(1+12*npsp:13*npsp) = pspheads(1:npsp)%nprojso(3)
107  list_int(1+13*npsp)         = test_paw
108 
109  call xmpi_bcast(list_int,master,comm,ierr)
110 
111  pspheads(1:npsp)%nproj(0)   = list_int(1        :   npsp)
112  pspheads(1:npsp)%nproj(1)   = list_int(1+   npsp: 2*npsp)
113  pspheads(1:npsp)%nproj(2)   = list_int(1+ 2*npsp: 3*npsp)
114  pspheads(1:npsp)%nproj(3)   = list_int(1+ 3*npsp: 4*npsp)
115  pspheads(1:npsp)%lmax       = list_int(1+ 4*npsp: 5*npsp)
116  pspheads(1:npsp)%xccc       = list_int(1+ 5*npsp: 6*npsp)
117  pspheads(1:npsp)%pspxc      = list_int(1+ 6*npsp: 7*npsp)
118  pspheads(1:npsp)%pspdat     = list_int(1+ 7*npsp: 8*npsp)
119  pspheads(1:npsp)%pspcod     = list_int(1+ 8*npsp: 9*npsp)
120  pspheads(1:npsp)%pspso      = list_int(1+ 9*npsp:10*npsp)
121  pspheads(1:npsp)%nprojso(1) = list_int(1+10*npsp:11*npsp)
122  pspheads(1:npsp)%nprojso(2) = list_int(1+11*npsp:12*npsp)
123  pspheads(1:npsp)%nprojso(3) = list_int(1+12*npsp:13*npsp)
124  test_paw                    = list_int(1+13*npsp)
125  ABI_DEALLOCATE(list_int)
126 
127 !Unbeliveable, this cannot be sent with the others, for woopy
128  ABI_ALLOCATE(list_int,(npsp))
129  list_int(1:npsp) = pspheads(1:npsp)%usewvl
130  call xmpi_bcast(list_int,master,comm,ierr)
131  pspheads(1:npsp)%usewvl     = list_int(1:npsp)
132  ABI_DEALLOCATE(list_int)
133 
134 
135 !Broadcast zionpsp and znuclpsp
136  ABI_ALLOCATE(list_dpr,(7*npsp))
137  list_dpr(1       :  npsp) = pspheads(1:npsp)%zionpsp
138  list_dpr(1+  npsp:2*npsp) = pspheads(1:npsp)%znuclpsp
139  list_dpr(1+2*npsp:3*npsp) = pspheads(1:npsp)%GTHradii(0)
140  list_dpr(1+3*npsp:4*npsp) = pspheads(1:npsp)%GTHradii(1)
141  list_dpr(1+4*npsp:5*npsp) = pspheads(1:npsp)%GTHradii(2)
142  list_dpr(1+5*npsp:6*npsp) = pspheads(1:npsp)%GTHradii(3)
143  list_dpr(1+6*npsp:7*npsp) = pspheads(1:npsp)%GTHradii(4)
144  
145  call xmpi_bcast(list_dpr,master,comm,ierr)
146 
147  pspheads(1:npsp)%zionpsp     = list_dpr(1       :  npsp)
148  pspheads(1:npsp)%znuclpsp    = list_dpr(1+  npsp:2*npsp)
149  pspheads(1:npsp)%GTHradii(0) = list_dpr(1+2*npsp:3*npsp)
150  pspheads(1:npsp)%GTHradii(1) = list_dpr(1+3*npsp:4*npsp)
151  pspheads(1:npsp)%GTHradii(2) = list_dpr(1+4*npsp:5*npsp)
152  pspheads(1:npsp)%GTHradii(3) = list_dpr(1+5*npsp:6*npsp)
153  pspheads(1:npsp)%GTHradii(4) = list_dpr(1+6*npsp:7*npsp)
154  ABI_DEALLOCATE(list_dpr)
155 
156 !Broadcast additional integers for PAW psps (testpaw was spread, previously)
157  if (test_paw==1) then
158    ABI_ALLOCATE(list_int,(6*npsp))
159    list_int(1       :  npsp)=pspheads(1:npsp)%pawheader%basis_size
160    list_int(1+  npsp:2*npsp)=pspheads(1:npsp)%pawheader%l_size
161    list_int(1+2*npsp:3*npsp)=pspheads(1:npsp)%pawheader%lmn_size
162    list_int(1+3*npsp:4*npsp)=pspheads(1:npsp)%pawheader%mesh_size
163    list_int(1+4*npsp:5*npsp)=pspheads(1:npsp)%pawheader%pawver
164    list_int(1+5*npsp:6*npsp)=pspheads(1:npsp)%pawheader%shape_type
165 
166    call xmpi_bcast(list_int,master,comm,ierr)
167 
168    pspheads(1:npsp)%pawheader%basis_size=list_int(1       :  npsp)
169    pspheads(1:npsp)%pawheader%l_size    =list_int(1+  npsp:2*npsp)
170    pspheads(1:npsp)%pawheader%lmn_size  =list_int(1+2*npsp:3*npsp)
171    pspheads(1:npsp)%pawheader%mesh_size =list_int(1+3*npsp:4*npsp)
172    pspheads(1:npsp)%pawheader%pawver    =list_int(1+4*npsp:5*npsp)
173    pspheads(1:npsp)%pawheader%shape_type=list_int(1+5*npsp:6*npsp)
174    ABI_DEALLOCATE(list_int)
175 
176 !  broadcast rpaw values
177    ABI_ALLOCATE(list_dpr,(2*npsp))
178 
179    list_dpr(1       :  npsp) = pspheads(1:npsp)%pawheader%rpaw
180    list_dpr(1+1*npsp:2*npsp) = pspheads(1:npsp)%pawheader%rshp
181 
182    call xmpi_bcast(list_dpr,master,comm,ierr)
183 
184    pspheads(1:npsp)%pawheader%rpaw = list_dpr(1       :  npsp)
185    pspheads(1:npsp)%pawheader%rshp = list_dpr(1+  npsp:2*npsp)
186 
187    ABI_DEALLOCATE(list_dpr)
188  end if
189 
190  call timab(48,2,tsec)
191 
192 #else
193 !Code to use unused dummy arguments
194  if(pspheads(1)%lmax == -10) pspheads(1)%lmax=-10
195  if(test_paw == -1) test_paw = -1
196 #endif
197 
198 end subroutine pspheads_comm