1!
2! Copyright (C) 2001-2020 Quantum ESPRESSO group
3! This file is distributed under the terms of the
4! GNU General Public License. See the file `License'
5! in the root directory of the present distribution,
6! or http://www.gnu.org/copyleft/gpl.txt .
7!
8!=----------------------------------------------------------------------------=!
9MODULE read_pseudo_mod
10!=----------------------------------------------------------------------------=!
11  !
12  !! read pseudopotential files and store the data on internal variables of the
13  !! program. Note that all processors read the same file!
14  !
15  USE io_files,     ONLY: pseudo_dir, pseudo_dir_cur, psfile, tmp_dir
16  USE ions_base,    ONLY: ntyp => nsp
17  !! global variables  required on input
18  !
19  USE atom,         ONLY: msh, rgrid
20  USE ions_base,    ONLY: zv
21  USE uspp_param,   ONLY: upf, nvb
22  USE uspp,         ONLY: okvan, nlcc_any
23  !! global variables modified on output
24  !
25  IMPLICIT NONE
26  SAVE
27  PRIVATE
28  !
29  PUBLIC :: readpp, check_order
30  !
31  CONTAINS
32  !
33  !-----------------------------------------------------------------------
34SUBROUTINE readpp ( input_dft, printout, ecutwfc_pp, ecutrho_pp )
35  !-----------------------------------------------------------------------
36  !
37  !! Reads PP files and puts the result into the "upf" structure of module uspp_param
38  !! Sets  DFT to input_dft if present, to the value read in PP files otherwise
39  !! Sets  number of valence electrons Zv, control variables okvan and nlcc_any,
40  !! compatibility variable nvb
41  !! Optionally returns cutoffs read from PP files into ecutwfc_pp, ecutrho_pp
42  !
43  USE kinds,        ONLY: DP
44  USE mp,           ONLY: mp_bcast, mp_sum
45  USE mp_images,    ONLY: intra_image_comm
46  USE io_global,    ONLY: stdout, ionode, ionode_id
47  USE pseudo_types, ONLY: pseudo_upf, deallocate_pseudo_upf
48  USE funct,        ONLY: enforce_input_dft, set_dft_from_name, &
49       get_iexch, get_icorr, get_igcx, get_igcc, get_inlc
50  use radial_grids, ONLY: deallocate_radial_grid, nullify_radial_grid
51  USE wrappers,     ONLY: md5_from_file, f_remove
52  USE read_upf_v1_module,   ONLY: read_upf_v1
53#if defined (__use_fox)
54  USE upf_module,   ONLY: read_upf_new
55  USE emend_upf_module, ONLY: make_emended_upf_copy
56#else
57  USE read_upf_new_module,  ONLY: read_upf_new
58#endif
59  USE upf_auxtools, ONLY: upf_get_pp_format, upf_check_atwfc_norm
60  USE upf_to_internal,  ONLY: add_upf_grid, set_upf_q
61  USE read_uspp_module, ONLY: readvan, readrrkj
62  USE m_gth,            ONLY: readgth
63  !
64  IMPLICIT NONE
65  !
66  CHARACTER(len=*), INTENT(INOUT) :: input_dft
67  LOGICAL, OPTIONAL, INTENT(IN) :: printout
68  REAL(DP), OPTIONAL, INTENT(OUT) :: ecutwfc_pp, ecutrho_pp
69  !
70  REAL(DP), parameter :: rcut = 10.d0
71  ! 2D Coulomb cutoff: modify this (at your own risks) if problems with cutoff
72  ! being smaller than pseudo rcut. original value=10.0
73  CHARACTER(len=512) :: file_pseudo ! file name complete with path
74  CHARACTER(len=512) :: file_fixed, msg
75  LOGICAL :: printout_ = .FALSE., exst, is_xml
76  INTEGER :: iunps, isupf, nt, nb, ir, ios
77  INTEGER :: iexch_, icorr_, igcx_, igcc_, inlc_
78  !
79  ! ... initializations, allocations, etc
80  !
81  iunps = 4
82  !
83  IF( ALLOCATED( upf ) ) THEN
84     DO nt = 1, SIZE( upf )
85        CALL deallocate_pseudo_upf( upf( nt ) )
86     END DO
87     DEALLOCATE( upf )
88  END IF
89  !
90  ALLOCATE ( upf( ntyp ) )
91  !
92  IF ( PRESENT(printout) ) THEN
93     printout_ = printout .AND. ionode
94  END IF
95  IF ( printout_) THEN
96     WRITE( stdout,"(//,3X,'Atomic Pseudopotentials Parameters',/, &
97                   &    3X,'----------------------------------' )" )
98  END IF
99  !
100  DO nt = 1, ntyp
101     !
102     ! try first pseudo_dir_cur if set: in case of restart from file,
103     ! this is where PP files should be located
104     !
105     ios = 1
106     IF ( pseudo_dir_cur /= ' ' ) THEN
107        file_pseudo  = TRIM (pseudo_dir_cur) // TRIM (psfile(nt))
108        INQUIRE(file = file_pseudo, EXIST = exst)
109        IF (exst) ios = 0
110        CALL mp_sum (ios,intra_image_comm)
111        IF ( ios /= 0 ) CALL infomsg &
112                     ('readpp', 'file '//TRIM(file_pseudo)//' not found')
113     END IF
114     !
115     ! file not found? no panic (yet): try the original location pseudo_dir
116     ! as set in input (it should already contain a slash at the end)
117     !
118     IF ( ios /= 0 ) THEN
119        file_pseudo = TRIM (pseudo_dir) // TRIM (psfile(nt))
120        INQUIRE ( file = file_pseudo, EXIST = exst)
121        IF (exst) ios = 0
122        CALL mp_sum (ios,intra_image_comm)
123        CALL errore('readpp', 'file '//TRIM(file_pseudo)//' not found',ABS(ios))
124     END IF
125     !
126     IF( printout_ ) THEN
127        WRITE( stdout, "(/,3X,'Reading pseudopotential for specie # ',I2, &
128                       & ' from file :',/,3X,A)") nt, TRIM(file_pseudo)
129     END IF
130     !
131     IF ( ionode ) THEN
132        isupf = 0
133        CALL  read_upf_new( file_pseudo, upf(nt), isupf )
134        !
135        !! start reading - check  first if files are readable as xml files,
136        !! then as UPF v.2, then as UPF v.1
137        !
138        IF (isupf ==-81 ) THEN
139#if defined (__use_fox)
140           !! error -81 may mean that file contains offending characters
141           !! fix and write file to tmp_dir
142           !! the underscore is added to distinguish this "fixed" file
143           !! from the original one, in case the latter is in tmp_dir
144           !
145           file_fixed = TRIM(tmp_dir)//TRIM(psfile(nt))//'_'
146           is_xml = make_emended_upf_copy( file_pseudo, file_fixed )
147           !
148           IF (is_xml) THEN
149              !
150              CALL  read_upf_new( file_fixed, upf(nt), isupf )
151              !! try again to read from the corrected file
152              WRITE ( msg, '(A)') 'Pseudo file '// trim(psfile(nt)) // ' has been fixed on the fly.' &
153            &    // new_line('a') // '     To avoid this message in the future, permanently fix ' &
154            &    // new_line('a') // '     your pseudo files following these instructions: ' &
155            &    // new_line('a') // '     https://gitlab.com/QEF/q-e/blob/master/upftools/how_to_fix_upf.md'
156             CALL infomsg('read_upf', trim(msg) )
157           ELSE
158              !
159              CALL  read_upf_v1 (file_pseudo, upf(nt), isupf )
160              !! try to read UPF v.1 file
161              IF ( isupf == 0 ) isupf = -1
162              !
163           END IF
164           !
165           ios = f_remove( file_fixed )
166#else
167           CALL  read_upf_v1 (file_pseudo, upf(nt), isupf )
168           !! try to read UPF v.1 file
169           IF ( isupf == 0 ) isupf = -1
170#endif
171        END IF
172        !
173     END IF
174     !
175     CALL mp_bcast (isupf,ionode_id,intra_image_comm)
176     !
177     IF (isupf == -2 .OR. isupf == -1 .OR. isupf == 0) THEN
178        !
179        CALL upf_bcast(upf(nt), ionode, ionode_id, intra_image_comm)
180        !! broadcast the pseudopotential to all processors
181        !
182        IF( printout_) THEN
183           IF ( isupf == 0 ) THEN
184              WRITE( stdout, "(3X,'file type is xml')")
185           ELSE
186              WRITE( stdout, "(3X,'file type is UPF v.',I1)") ABS(isupf)
187           END IF
188        END IF
189        !
190     ELSE
191        !
192        OPEN ( UNIT = iunps, FILE = file_pseudo, STATUS = 'old', FORM = 'formatted' )
193        !
194        !     The type of the pseudopotential is determined by the file name:
195        !    *.xml or *.XML  UPF format with schema              pp_format=0
196        !    *.upf or *.UPF  UPF format                          pp_format=1
197        !    *.vdb or *.van  Vanderbilt US pseudopotential code  pp_format=2
198        !    *.gth           Goedecker-Teter-Hutter NC pseudo    pp_format=3
199        !    *.RRKJ3         Andrea's   US new code              pp_format=4
200        !    none of the above: PWSCF norm-conserving format     pp_format=5
201        !
202        IF ( upf_get_pp_format( psfile(nt) ) == 2  ) THEN
203           !
204           IF( printout_ ) &
205              WRITE( stdout, "(3X,'file type is Vanderbilt US PP')")
206           CALL readvan (iunps, nt, upf(nt))
207           !
208        ELSE IF ( upf_get_pp_format( psfile(nt) ) == 3 ) THEN
209           !
210           IF( printout_ ) &
211              WRITE( stdout, "(3X,'file type is GTH (analytical)')")
212           CALL readgth (iunps, nt, upf(nt))
213           !
214        ELSE IF ( upf_get_pp_format( psfile(nt) ) == 4 ) THEN
215           !
216           IF( printout_ ) &
217              WRITE( stdout, "(3X,'file type is RRKJ3')")
218           CALL readrrkj (iunps, nt, upf(nt))
219           !
220        ELSE IF ( upf_get_pp_format( psfile(nt) ) == 5 ) THEN
221           !
222           IF( printout_ ) &
223              WRITE( stdout, "(3X,'file type is old PWscf NC format')")
224           CALL read_ncpp (iunps, nt, upf(nt))
225           !
226        ELSE
227           !
228           CALL errore('readpp', 'file '//TRIM(file_pseudo)//' not readable',1)
229           !
230        ENDIF
231        !
232        ! end of reading
233        !
234        CLOSE (iunps)
235        !
236     ENDIF
237     !
238     ! reconstruct Q(r) if needed
239     !
240     CALL set_upf_q (upf(nt))
241     !
242     ! Calculate MD5 checksum for this pseudopotential
243     !
244     CALL md5_from_file(file_pseudo, upf(nt)%md5_cksum)
245     !
246  END DO
247  !
248  ! end of PP reading - now set up more variables
249  !
250  ! radial grids -
251  !
252  IF( ALLOCATED( rgrid ) ) THEN
253     DO nt = 1, SIZE( rgrid )
254        CALL deallocate_radial_grid( rgrid( nt ) )
255        CALL nullify_radial_grid( rgrid( nt ) )
256     END DO
257     DEALLOCATE( rgrid )
258     if(allocated(msh)) DEALLOCATE( msh )
259  END IF
260  ALLOCATE( rgrid( ntyp ), msh( ntyp ) )
261  !
262  nvb = 0
263  DO nt = 1, ntyp
264     !
265     CALL nullify_radial_grid( rgrid( nt ) )
266     CALL add_upf_grid (upf(nt), rgrid(nt))
267     !
268     ! the radial grid is defined up to r(mesh) but we introduce
269     ! an auxiliary variable msh to limit the grid up to rcut=10 a.u.
270     ! This is used to cut off the numerical noise arising from the
271     ! large-r tail in cases like the integration of V_loc-Z/r
272     !
273     DO ir = 1, rgrid(nt)%mesh
274        IF (rgrid(nt)%r(ir) > rcut) THEN
275           msh (nt) = ir
276           GOTO 5
277        END IF
278     END DO
279     msh (nt) = rgrid(nt)%mesh
2805    msh (nt) = 2 * ( (msh (nt) + 1) / 2) - 1
281     !
282     ! msh is forced to be odd for simpson integration (maybe obsolete?)
283     !
284     ! ... Zv = valence charge of the (pseudo-)atom, read from PP files,
285     ! ... is set equal to Zp = pseudo-charge of the pseudopotential
286     !
287     zv(nt) = upf(nt)%zp
288     !
289     ! ... count US species (obsolete?)
290     !
291     IF (upf(nt)%tvanp) nvb=nvb+1
292     !
293     ! check for zero atomic wfc,
294     ! check that (occupied) atomic wfc are properly normalized
295     !
296     CALL upf_check_atwfc_norm(upf(nt),psfile(nt))
297     !
298  END DO
299  !
300  ! ... set DFT value
301  !
302  IF (input_dft /='none') CALL enforce_input_dft (input_dft)
303  !
304  DO nt = 1, ntyp
305     !
306     CALL set_dft_from_name( upf(nt)%dft )
307     !
308     ! ... Check for DFT consistency - ignored if dft enforced from input
309     !
310     IF (nt == 1) THEN
311        iexch_ = get_iexch()
312        icorr_ = get_icorr()
313        igcx_  = get_igcx()
314        igcc_  = get_igcc()
315        inlc_  = get_inlc()
316     ELSE
317        IF ( iexch_ /= get_iexch() .OR. icorr_ /= get_icorr() .OR. &
318             igcx_  /= get_igcx()  .OR. igcc_  /= get_igcc()  .OR.  &
319             inlc_  /= get_inlc() ) THEN
320           CALL errore( 'readpp','inconsistent DFT read from PP files', nt)
321        END IF
322     END IF
323     !
324  END DO
325  !
326  ! more initializations
327  !
328  okvan = ( nvb > 0 )
329  nlcc_any = ANY ( upf(1:ntyp)%nlcc )
330  !
331  ! return cutoff read from PP file, if required
332  !
333  IF ( PRESENT(ecutwfc_pp) ) THEN
334     ecutwfc_pp = MAXVAL ( upf(1:ntyp)%ecutwfc )
335  END IF
336  IF ( PRESENT(ecutrho_pp) ) THEN
337     ecutrho_pp = MAXVAL ( upf(1:ntyp)%ecutrho )
338  END IF
339  !
340END SUBROUTINE readpp
341!
342SUBROUTINE check_order
343   ! CP-specific check
344   IF ( ANY(upf(1:ntyp)%tpawp) ) CALL errore ('readpp','PAW not implemented',1)
345END SUBROUTINE check_order
346!
347! Copyright (C) 2020 Quantum ESPRESSO group
348! This file is distributed under the terms of the
349! GNU General Public License. See the file `License'
350! in the root directory of the present distribution,
351! or http://www.gnu.org/copyleft/gpl.txt .
352!
353!------------------------------------------------+
354SUBROUTINE upf_bcast(upf, ionode, ionode_id, comm)
355  !---------------------------------------------+
356  !
357  !! Broadcast the "upf" structure, read on processor "ionode_id",
358  !! to all other processors in the communicator "comm".
359  !
360  USE kinds,        ONLY: DP
361  USE pseudo_types, ONLY: pseudo_upf
362  USE mp,           ONLY: mp_bcast
363  !
364  IMPLICIT NONE
365  !
366  TYPE(pseudo_upf),INTENT(INOUT) :: upf
367  !! pseudo_upf type structure storing the pseudo data
368  LOGICAL, INTENT(in) :: ionode
369  !! true if we are on the processor that broadcasts
370  !! upf is allocated if (ionode), must be allocated otherwise
371  INTEGER, INTENT(in) :: ionode_id
372  !! ID of the processor that broadcasts
373  INTEGER, INTENT(in) :: comm
374  !! MPI communicator
375  !
376  CALL mp_bcast (upf%nv, ionode_id, comm )
377  CALL mp_bcast (upf%generated, ionode_id, comm )
378  CALL mp_bcast (upf%author, ionode_id, comm )
379  CALL mp_bcast (upf%date, ionode_id, comm )
380  CALL mp_bcast (upf%comment, ionode_id, comm )
381  CALL mp_bcast (upf%psd, ionode_id, comm )
382  CALL mp_bcast (upf%typ, ionode_id, comm )
383  CALL mp_bcast (upf%rel, ionode_id, comm )
384  CALL mp_bcast (upf%tvanp, ionode_id, comm )
385  CALL mp_bcast (upf%tpawp, ionode_id, comm )
386  CALL mp_bcast (upf%tcoulombp, ionode_id, comm )
387  CALL mp_bcast (upf%is_gth, ionode_id, comm )
388  CALL mp_bcast (upf%is_multiproj, ionode_id, comm )
389  CALL mp_bcast (upf%has_so, ionode_id, comm )
390  CALL mp_bcast (upf%has_wfc, ionode_id, comm )
391  CALL mp_bcast (upf%has_gipaw, ionode_id, comm )
392  CALL mp_bcast (upf%paw_as_gipaw, ionode_id, comm )
393  CALL mp_bcast (upf%nlcc, ionode_id, comm )
394  CALL mp_bcast (upf%dft, ionode_id, comm )
395  CALL mp_bcast (upf%zp, ionode_id, comm )
396  CALL mp_bcast (upf%etotps, ionode_id, comm )
397  CALL mp_bcast (upf%ecutwfc, ionode_id, comm )
398  CALL mp_bcast (upf%ecutrho, ionode_id, comm )
399  CALL mp_bcast (upf%lmax, ionode_id, comm )
400  CALL mp_bcast (upf%lmax_rho, ionode_id, comm )
401  CALL mp_bcast (upf%lloc, ionode_id, comm )
402  CALL mp_bcast (upf%mesh, ionode_id, comm )
403  CALL mp_bcast (upf%nwfc, ionode_id, comm )
404  CALL mp_bcast (upf%nbeta, ionode_id, comm )
405  CALL mp_bcast (upf%dx, ionode_id, comm )
406  CALL mp_bcast (upf%xmin, ionode_id, comm )
407  CALL mp_bcast (upf%rmax, ionode_id, comm )
408  CALL mp_bcast (upf%zmesh, ionode_id, comm )
409  !
410  IF ( .NOT. ionode) ALLOCATE( upf%r( upf%mesh ), upf%rab( upf%mesh ) )
411  CALL mp_bcast (upf%r,   ionode_id, comm )
412  CALL mp_bcast (upf%rab, ionode_id, comm )
413  !
414  IF ( .NOT. ionode) ALLOCATE( upf%rho_atc(upf%mesh) )
415  CALL mp_bcast (upf%rho_atc, ionode_id, comm )
416  !
417  IF(.not. upf%tcoulombp) THEN
418     IF ( .NOT. ionode) ALLOCATE( upf%vloc(upf%mesh) )
419     CALL mp_bcast (upf%vloc, ionode_id, comm )
420  ENDIF
421  !
422  IF ( .not. ionode) THEN
423     IF ( upf%nbeta == 0) THEN
424        upf%nqf = 0
425        upf%nqlc= 0
426        upf%qqq_eps= -1._dp
427        upf%kkbeta = 0
428        ALLOCATE( upf%kbeta(1),     &
429             upf%lll(1),           &
430             upf%beta(upf%mesh,1), &
431             upf%dion(1,1),        &
432             upf%rcut(1),          &
433             upf%rcutus(1),        &
434             upf%els_beta(1) )
435     ELSE
436        ALLOCATE( upf%kbeta(upf%nbeta),     &
437             upf%lll(upf%nbeta),            &
438             upf%beta(upf%mesh, upf%nbeta), &
439             upf%dion(upf%nbeta, upf%nbeta),&
440             upf%rcut(upf%nbeta),           &
441             upf%rcutus(upf%nbeta),         &
442             upf%els_beta(upf%nbeta) )
443     END IF
444  END IF
445  !
446  CALL mp_bcast (upf%beta, ionode_id, comm )
447  CALL mp_bcast (upf%kbeta, ionode_id, comm )
448  CALL mp_bcast (upf%els_beta, ionode_id, comm )
449  CALL mp_bcast (upf%lll, ionode_id, comm )
450  CALL mp_bcast (upf%rcut, ionode_id, comm )
451  CALL mp_bcast (upf%rcutus, ionode_id, comm )
452  CALL mp_bcast (upf%dion, ionode_id, comm )
453
454  IF(upf%tvanp .or. upf%tpawp) THEN
455     CALL mp_bcast (upf%q_with_l, ionode_id, comm )
456     CALL mp_bcast (upf%nqf, ionode_id, comm )
457     CALL mp_bcast (upf%nqlc, ionode_id, comm )
458     IF (upf%tpawp) THEN
459        IF ( .not. ionode) ALLOCATE &
460             ( upf%paw%augmom(upf%nbeta,upf%nbeta, 0:2*upf%lmax) )
461        CALL mp_bcast (upf%paw%augshape, ionode_id, comm )
462        CALL mp_bcast (upf%paw%raug, ionode_id, comm )
463        CALL mp_bcast (upf%paw%iraug, ionode_id, comm )
464        CALL mp_bcast (upf%paw%lmax_aug, ionode_id, comm )
465        CALL mp_bcast (upf%paw%augmom, ionode_id, comm )
466     END IF
467     CALL mp_bcast (upf%qqq_eps, ionode_id, comm )
468     IF ( .not. ionode) THEN
469        IF ( upf%nbeta == 0 ) THEN
470           ALLOCATE(upf%rinner(1), &
471             upf%qqq(1,1),         &
472             upf%qfunc(upf%mesh,1),&
473             upf%qfcoef(1,1,1,1) )
474             IF ( upf%q_with_l ) &
475                ALLOCATE( upf%qfuncl ( upf%mesh, 1, 1 ) )
476        ELSE
477           ALLOCATE( upf%qqq   ( upf%nbeta, upf%nbeta ) )
478           IF ( upf%q_with_l ) THEN
479              ALLOCATE( upf%qfuncl ( upf%mesh, upf%nbeta*(upf%nbeta+1)/2, 0:2*upf%lmax ) )
480           ELSE
481              ALLOCATE( upf%qfunc (upf%mesh, upf%nbeta*(upf%nbeta+1)/2) )
482           ENDIF
483           ALLOCATE( upf%rinner( upf%nqlc ) )
484           IF(upf%nqf <= 0) THEN
485              ALLOCATE( upf%qfcoef(1,1,1,1) )
486           ELSE
487              ALLOCATE( upf%qfcoef( upf%nqf, upf%nqlc, &
488                   upf%nbeta, upf%nbeta ) )
489           END IF
490        END IF
491     ENDIF
492     CALL mp_bcast (upf%qqq   , ionode_id, comm )
493     CALL mp_bcast (upf%rinner, ionode_id, comm )
494     CALL mp_bcast (upf%qfcoef, ionode_id, comm )
495     IF (upf%q_with_l) THEN
496        CALL mp_bcast (upf%qfuncl, ionode_id, comm )
497     ELSE
498        CALL mp_bcast (upf%qfunc , ionode_id, comm )
499     END IF
500     !
501  END IF
502  upf%kkbeta = MAXVAL(upf%kbeta(1:upf%nbeta))
503  IF(upf%tpawp) upf%kkbeta = MAX(upf%kkbeta, upf%paw%iraug)
504
505  IF ( .not. ionode ) THEN
506     ALLOCATE( upf%chi(upf%mesh,upf%nwfc) )
507     ALLOCATE( upf%els(upf%nwfc), &
508          upf%oc(upf%nwfc), &
509          upf%lchi(upf%nwfc), &
510          upf%nchi(upf%nwfc), &
511          upf%rcut_chi(upf%nwfc), &
512          upf%rcutus_chi(upf%nwfc), &
513          upf%epseu(upf%nwfc) )
514  END IF
515  CALL mp_bcast (upf%chi,ionode_id, comm )
516  CALL mp_bcast (upf%els, ionode_id, comm )
517  CALL mp_bcast (upf%oc,ionode_id, comm )
518  CALL mp_bcast (upf%lchi,ionode_id, comm )
519  CALL mp_bcast (upf%nchi,ionode_id, comm )
520  CALL mp_bcast (upf%rcut_chi,ionode_id, comm )
521  CALL mp_bcast (upf%rcutus_chi,ionode_id, comm )
522  CALL mp_bcast (upf%epseu,ionode_id, comm )
523  !
524  IF(upf%has_wfc) THEN
525     IF ( .not. ionode) THEN
526        ALLOCATE( upf%aewfc(upf%mesh, upf%nbeta) )
527        ALLOCATE( upf%pswfc(upf%mesh, upf%nbeta) )
528        IF (upf%has_so .and. upf%tpawp) ALLOCATE &
529             ( upf%paw%aewfc_rel(upf%mesh, upf%nbeta) )
530     END IF
531     IF (upf%has_so .and. upf%tpawp) CALL mp_bcast &
532          (upf%paw%aewfc_rel,ionode_id,comm )
533     CALL mp_bcast &
534          (upf%aewfc,ionode_id,comm )
535     CALL mp_bcast &
536          (upf%pswfc,ionode_id,comm )
537  END IF
538  !
539  IF ( .not. ionode) ALLOCATE( upf%rho_at(upf%mesh) )
540  CALL mp_bcast (upf%rho_at,ionode_id,comm )
541
542  IF (upf%has_so) THEN
543     IF ( .NOT. ionode) THEN
544        ALLOCATE (upf%nn(upf%nwfc))
545        ALLOCATE (upf%jchi(upf%nwfc))
546        ALLOCATE(upf%jjj(upf%nbeta))
547     END IF
548     CALL mp_bcast (upf%nn,ionode_id,comm )
549     CALL mp_bcast (upf%jchi,ionode_id,comm )
550     CALL mp_bcast (upf%jjj,ionode_id,comm )
551  END IF
552
553  IF (upf%tpawp) THEN
554     CALL mp_bcast (upf%paw_data_format,ionode_id,comm )
555     CALL mp_bcast (upf%paw%core_energy,ionode_id,comm )
556     IF ( .not. ionode ) THEN
557        ALLOCATE( upf%paw%oc(upf%nbeta) )
558        ALLOCATE( upf%paw%ae_rho_atc(upf%mesh) )
559        ALLOCATE( upf%paw%ae_vloc(upf%mesh) )
560        ALLOCATE( upf%paw%pfunc(upf%mesh, upf%nbeta,upf%nbeta) )
561        ALLOCATE(upf%paw%ptfunc(upf%mesh, upf%nbeta,upf%nbeta) )
562        IF (upf%has_so) &
563             ALLOCATE(upf%paw%pfunc_rel(upf%mesh, upf%nbeta,upf%nbeta) )
564     END IF
565     CALL mp_bcast (upf%paw%oc,ionode_id,comm )
566     CALL mp_bcast (upf%paw%ae_rho_atc,ionode_id,comm )
567     CALL mp_bcast (upf%paw%ae_vloc,ionode_id,comm )
568     CALL mp_bcast (upf%paw%pfunc,ionode_id,comm )
569     CALL mp_bcast (upf%paw%ptfunc,ionode_id,comm )
570     IF (upf%has_so) &
571          CALL mp_bcast (upf%paw%pfunc_rel,ionode_id,comm )
572  END IF
573
574  IF (upf%has_gipaw) THEN
575     CALL mp_bcast (upf%gipaw_data_format,ionode_id,comm )
576     CALL mp_bcast (upf%gipaw_ncore_orbitals,ionode_id,comm )
577     IF ( .not. ionode) THEN
578        ALLOCATE ( upf%gipaw_core_orbital_n(upf%gipaw_ncore_orbitals) )
579        ALLOCATE ( upf%gipaw_core_orbital_el(upf%gipaw_ncore_orbitals) )
580        ALLOCATE ( upf%gipaw_core_orbital_l(upf%gipaw_ncore_orbitals) )
581        ALLOCATE ( upf%gipaw_core_orbital(upf%mesh,upf%gipaw_ncore_orbitals) )
582     END IF
583     CALL mp_bcast (upf%gipaw_core_orbital_n ,ionode_id,comm )
584     CALL mp_bcast (upf%gipaw_core_orbital_el,ionode_id,comm )
585     CALL mp_bcast (upf%gipaw_core_orbital_l ,ionode_id,comm )
586     CALL mp_bcast (upf%gipaw_core_orbital   ,ionode_id,comm )
587     CALL mp_bcast (upf%gipaw_wfs_nchannels  ,ionode_id,comm )
588     IF ( .not. ionode) THEN
589        ALLOCATE ( upf%gipaw_wfs_el(upf%gipaw_wfs_nchannels) )
590        ALLOCATE ( upf%gipaw_wfs_ll(upf%gipaw_wfs_nchannels) )
591        ALLOCATE ( upf%gipaw_wfs_rcut(upf%gipaw_wfs_nchannels) )
592        ALLOCATE ( upf%gipaw_wfs_rcutus(upf%gipaw_wfs_nchannels) )
593        ALLOCATE ( upf%gipaw_wfs_ae(upf%mesh,upf%gipaw_wfs_nchannels) )
594        ALLOCATE ( upf%gipaw_wfs_ps(upf%mesh,upf%gipaw_wfs_nchannels) )
595        ALLOCATE ( upf%gipaw_vlocal_ae(upf%mesh) )
596        ALLOCATE ( upf%gipaw_vlocal_ps(upf%mesh) )
597     END IF
598     CALL mp_bcast (upf%gipaw_wfs_el, ionode_id,comm )
599     CALL mp_bcast (upf%gipaw_wfs_ll, ionode_id,comm )
600     CALL mp_bcast (upf%gipaw_wfs_rcut  , ionode_id,comm )
601     CALL mp_bcast (upf%gipaw_wfs_rcutus, ionode_id,comm )
602     CALL mp_bcast (upf%gipaw_wfs_ae,ionode_id,comm )
603     CALL mp_bcast (upf%gipaw_wfs_ps,ionode_id,comm )
604     CALL mp_bcast (upf%gipaw_vlocal_ae,ionode_id,comm )
605     CALL mp_bcast (upf%gipaw_vlocal_ps,ionode_id,comm )
606     !
607  END IF
608  !
609END SUBROUTINE upf_bcast
610
611END MODULE read_pseudo_mod
612
613