TABLE OF CONTENTS


ABINIT/psp_dump_outputs [ Functions ]

[ Top ] [ Functions ]

NAME

 psp_dump_outputs

FUNCTION

 (To be described ...)

COPYRIGHT

 Copyright (C) 2017-2018 ABINIT group (YP)
 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

 (to be filled)

OUTPUT

 (to be filled)

SIDE EFFECTS

 (to be filled)

PARENTS

      pspatm

CHILDREN

SOURCE

 32 !      call psp_dump_outputs(psps%lmnmax,psps%mpssoang,psps%lnmax, &
 33 !&      psps%mqgrid_ff,psps%n1xccc,mmax,epsatm,qchrg,xcccrc,nctab, &
 34 !&      indlmn,nproj,ekb,ffspl,vlspl,xccc1d)
 35 
 36 #if defined HAVE_CONFIG_H
 37 #include "config.h"
 38 #endif
 39 
 40 #include "abi_common.h"
 41 
 42 subroutine psp_dump_outputs(pfx,pspcod,lmnmax,lnmax,mpssoang, &
 43 &      mqgrid,n1xccc,mmax,maxrad,epsatm,qchrg,xcccrc,nctab, &
 44 &      indlmn,nproj,ekb,ffspl,vlspl,xccc1d)
 45 
 46  use defs_basis
 47  use defs_datatypes, only : nctab_t
 48  use m_errors
 49 
 50 !This section has been created automatically by the script Abilint (TD).
 51 !Do not modify the following lines by hand.
 52 #undef ABI_FUNC
 53 #define ABI_FUNC 'psp_dump_outputs'
 54 !End of the abilint section
 55 
 56  implicit none
 57 
 58 !Arguments ------------------------------------
 59 !scalars
 60  character(len=*), intent(in) :: pfx
 61  integer,intent(in) :: pspcod,lmnmax,lnmax,mpssoang,mqgrid,n1xccc
 62  integer,intent(in) :: mmax
 63  real(dp),intent(in) :: maxrad,epsatm,qchrg,xcccrc
 64  type(nctab_t),intent(in) :: nctab
 65 !arrays
 66  integer,intent(in) :: indlmn(6,lmnmax),nproj(mpssoang)
 67  real(dp),intent(in) :: ekb(lnmax),ffspl(mqgrid,2,lnmax),vlspl(mqgrid,2)
 68  real(dp),intent(in) :: xccc1d(n1xccc,6)
 69 
 70 !Local variables ------------------------------
 71 !scalars
 72  integer, parameter :: dump = 64
 73  integer :: ierr, i, j ,k
 74  character(len=500) :: msg
 75 
 76  ! *********************************************************************
 77 
 78  open(unit=dump, file=trim(pfx)//"_psp_info.yaml", status='REPLACE', err=10, iostat=ierr)
 79 
 80  write(dump,'(3a)') "%YAML 1.2", ch10, "---"
 81 
 82  write(dump, '(2a)') ch10, "# Pseudopotential info"
 83  write(dump, '(a,1x,i8)') "pspcod:", pspcod
 84 
 85  write(dump, '(2a)') ch10, "# Array dimensions"
 86  write(dump, '(a)') "dims:"
 87  write(dump, '(4x,a,1x,i8)') "lmnmax:", lmnmax
 88  write(dump, '(4x,a,1x,i8)') "lnmax:", lnmax
 89  write(dump, '(4x,a,1x,i8)') "mpssoang:", mpssoang
 90  write(dump, '(4x,a,1x,i8)') "mqgrid:", mqgrid
 91  write(dump, '(4x,a,1x,i8)') "n1xccc:", n1xccc
 92  write(dump, '(4x,a,1x,i8)') "mmax:", mmax
 93 
 94  write(dump, '(2a)') ch10, "# Quantities"
 95  write(dump, '(a,1x,e12.5)') "maxrad:", maxrad
 96  write(dump, '(a,1x,e12.5)') "epsatm:", epsatm
 97  write(dump, '(a,1x,e12.5)') "qchrg:", qchrg
 98  write(dump, '(a,1x,e12.5)') "xcccrc:", xcccrc
 99 
100  write(dump, '(2a)') ch10, "# Structure: nctab"
101  write(dump, '(a)') "nctab:"
102  write(dump,'(4x,a,":",1x,i4)') "mqgrid_vl", nctab%mqgrid_vl
103  write(dump,'(4x,a,":",1x,l4)') "has_tvale", nctab%has_tvale
104  write(dump,'(4x,a,":",1x,l4)') "has_tcore", nctab%has_tcore
105  write(dump,'(4x,a,":",1x,e12.5)') "dncdq0", nctab%dncdq0
106  write(dump,'(4x,a,":",1x,e12.5)') "d2ncdq0", nctab%d2ncdq0
107  write(dump,'(4x,a,":",1x,e12.5)') "dnvdq0", nctab%dnvdq0
108 
109  if ( nctab%has_tvale ) then
110    write(dump, '(2a)') ch10, "# Array: nctab_tvalespl(mqgrid_vl,2)"
111    write(dump, '(a)') "nctab_tvalespl:"
112    do j=1,2
113      do i=1,nctab%mqgrid_vl
114        if ( i == 1 ) then
115          write(dump,'(4x,a,1x,e12.5)') "- -", nctab%tvalespl(i,j)
116        else
117          write(dump,'(4x,a,1x,e12.5)') "  -", nctab%tvalespl(i,j)
118        end if
119      end do
120    end do
121  end if
122 
123  if ( nctab%has_tcore ) then
124    write(dump, '(2a)') ch10, "# Array: nctab_tcorespl(mqgrid_vl,2)"
125    write(dump, '(a)') "nctab_tcorespl:"
126    do j=1,2
127      do i=1,nctab%mqgrid_vl
128        if ( i == 1 ) then
129          write(dump,'(4x,a,1x,e12.5)') "- -", nctab%tcorespl(i,j)
130        else
131          write(dump,'(4x,a,1x,e12.5)') "  -", nctab%tcorespl(i,j)
132        end if
133      end do
134    end do
135  end if
136 
137  write(dump, '(2a)') ch10, "# Array: integer indlmn(6,lmnmax)"
138  write(dump, '(a)') "indlmn:"
139  do i=1,lmnmax
140    write(dump,'(4x,a,i4,5(",",i4),a)') "- [", indlmn(:,i), "]"
141  end do
142 
143  write(dump, '(2a)') ch10, "# Array: integer nproj(mpssoang)"
144  write(dump, '(a)') "nproj:"
145  do i=1,mpssoang
146    write(dump,'(4x,"-",1x,i4)') nproj(i)
147  end do
148 
149  write(dump, '(2a)') ch10, "# Array: double ekb(lnmax)"
150  write(dump, '(a)') "ekb:"
151  do i=1,lnmax
152    write(dump,'(4x,"-",1x,e12.5)') ekb(i)
153  end do
154 
155  write(dump, '(2a)') ch10, "# Array: ffspl(mqgrid,2,lnmax)"
156  write(dump, '(a)') "ffspl:"
157  do k=1,lnmax
158    do j=1,2
159      do i=1,mqgrid
160        if ( (i == 1) .and. (j == 1) ) then
161          write(dump,'(4x,a,1x,e12.5)') "- - -", ffspl(i,j,k)
162        else if ( i == 1 ) then
163          write(dump,'(4x,a,1x,e12.5)') "  - -", ffspl(i,j,k)
164        else
165          write(dump,'(4x,a,1x,e12.5)') "    -", ffspl(i,j,k)
166        end if
167      end do
168    end do
169  end do
170 
171  write(dump, '(2a)') ch10, "# Array: vlspl(mqgrid,2)"
172  write(dump, '(a)') "vlspl:"
173  do j=1,2
174    do i=1,mqgrid
175      if ( i == 1 ) then
176        write(dump,'(4x,a,1x,e12.5)') "- -", vlspl(i,j)
177      else
178        write(dump,'(4x,a,1x,e12.5)') "  -", vlspl(i,j)
179      end if
180    end do
181  end do
182 
183  write(dump, '(2a)') ch10, "# Array: xccc1d(n1xccc,6)"
184  write(dump, '(a)') "xccc1d:"
185  do j=1,6
186    do i=1,mqgrid
187      if ( i == 1 ) then
188        write(dump,'(4x,a,1x,e12.5)') "- -", xccc1d(i,j)
189      else
190        write(dump,'(4x,a,1x,e12.5)') "  -", xccc1d(i,j)
191      end if
192    end do
193  end do
194 
195  write (dump,'(2a)') ch10, "..."
196 
197  close(dump)
198 
199  return
200 
201  10 continue
202 
203  if ( ierr /= 0 ) then
204    write(msg,'(a,a,a,i8)') "Error writing pseudopotential information", &
205 &   ch10, "IOSTAT=", ierr
206    MSG_WARNING(msg)
207  end if
208 
209 end subroutine psp_dump_outputs