3
subroutine pseudoXML( ray, npotd, npotu, zion, zratio )
15
include 'pseudowave.h'
20
double precision :: zion, zratio
22
character*4 :: polattrib, relattrib, coreattrib
23
character*10 :: ray(6)
24
character*30 xcfuntype, xcfunparam
25
character*30 gridtype, gridunits, gridscale, gridstep, gridnpoint
26
character*30 slpotunits, slpotformat, slpotndown, slpotnup
27
character*30 slwfnunits, slwfnformat, slwfndown, slwfnup
29
double precision, allocatable :: chval(:)
32
! These determine the format for ASCII files
34
character(len=*), parameter ::
35
. fmt_int = "(tr1,i2)" ,
36
. fmt_nam = "(tr1,a2,tr1,a2,tr1,a3,tr1,a4)",
37
. fmt_met = "(tr1,6a10,/,tr1,a70)" ,
38
. fmt_pot= "(tr1,2i3,i5,3g20.12)" ,
39
. fmt_rad = "(4(g20.12))" ,
43
! Digest and dump the information about the exchange and correlation functional
48
xcfunparam = 'Ceperley-Alder'
52
xcfunparam = 'Perdew-Wang-92'
60
xcfunparam = 'Hedin-Lundqvist'
64
xcfunparam = 'Gunnarson-Lundqvist'
68
xcfunparam = 'von Barth-Hedin'
72
xcfunparam = 'Perdew-Burke-Ernzerhof'
76
xcfunparam = 'RPBE - Hammer et al'
80
xcfunparam = 'revPBE Zhang+Yang'
84
xcfunparam = 'Becke-Lee-Yang-Parr'
88
xcfunparam = 'Wu-Cohen'
92
xcfunparam = 'Perdew-Burke-Ernzerhof-solid'
96
xcfunparam = 'Armiento-Mattsson-05'
100
xcfunparam = 'Perdew-Wang-91'
104
xcfunparam = 'Capelle-et-al--PBE-Js-Jr-LO'
108
xcfunparam = 'Capelle-et-al--PBE-Js-Jr-HEG'
112
xcfunparam = 'Capelle-et-al--PBE-Gc-Gx-LO'
116
xcfunparam = 'Capelle-et-al--PBE-Gc-Gx-HEG'
120
xcfunparam = 'Dion-et-al--DRSLL'
124
xcfunparam = 'Lee-et-al--LMKLL'
128
xcfunparam = 'Klimes-et-al--KBM'
132
xcfunparam = 'Cooper--C09'
136
xcfunparam = 'Berland-Hyldgaard--BH'
140
xcfunparam = 'Vydrov-VanVoorhis--VV'
144
! Digest and dump the information about the pseudopotential flavor
161
! Digest and dump the information about the non-linear core corrections
181
! Digest and dump the information about the grid
186
gridnpoint = str(nr-1)
188
! Digest and dump the information about the semilocal components
189
slpotunits = 'rydberg'
191
slpotndown = str(npotd)
192
slpotnup = str(npotu)
194
! Digest and dump the information about the pseudowave functions
195
slwfnunits = 'electrons/bohr^(-3/2)'
196
slwfnformat = 'u_n,l (r) = 1/r R_n,l (r)'
197
slwfndown = str(npotd)
200
! Allocate and define the valence charge density
201
allocate(chval(1:nr))
204
chval(ip) = zratio * (cdd(ip)+cdu(ip))
208
! ---------------------------------------------------------------------
210
call xml_OpenFile("VPSXML",xf, pretty_print=.true.)
212
! call xml_AddXMLDeclaration(xf,"UTF-8")
214
call xml_NewElement(xf,"pseudo")
216
call xml_NewElement(xf,"header")
217
call my_add_attribute(xf,"symbol",nameat)
218
call my_add_attribute(xf,"atomic-number",str(znuc))
219
call my_add_attribute(xf,"zval",str(zion))
220
call my_add_attribute(xf,"creator",ray(1))
221
call my_add_attribute(xf,"date",ray(2))
222
call my_add_attribute(xf,"flavor",ray(3)//ray(4))
223
call my_add_attribute(xf,"relativistic",relattrib)
224
call my_add_attribute(xf,"polarized",polattrib)
225
call my_add_attribute(xf,"core-corrections",coreattrib)
226
call my_add_attribute(xf,"xc-functional-type",xcfuntype)
227
call my_add_attribute(xf,"xc-functional-parametrization",
229
call xml_EndElement(xf,"header")
231
call xml_NewElement(xf,"grid")
232
call my_add_attribute(xf,"type",gridtype)
233
call my_add_attribute(xf,"units",gridunits)
234
call my_add_attribute(xf,"scale",gridscale)
235
call my_add_attribute(xf,"step",gridstep)
236
call my_add_attribute(xf,"npts",gridnpoint)
237
call xml_EndElement(xf,"grid")
239
call xml_NewElement(xf,"semilocal")
240
call my_add_attribute(xf,"units",slpotunits)
241
call my_add_attribute(xf,"format",slpotformat)
242
call my_add_attribute(xf,"npots-down",slpotndown)
243
call my_add_attribute(xf,"npots-up",slpotnup)
245
! Down pseudopotentials follows
247
vpsd: do ivps = 1, lmax
248
if (indd(ivps) .eq. 0) cycle
249
call xml_NewElement(xf,"vps")
250
call my_add_attribute(xf,"principal-n",str(no(indd(ivps))))
251
call my_add_attribute(xf,"l",il(ivps))
252
call my_add_attribute(xf,"cutoff",str(rc(ivps)))
253
call my_add_attribute(xf,"occupation",str(zo(indd(ivps))))
254
! JMS/AG 2009/04/02: to be checked which one works better
255
! call my_add_attribute(xf,"spin",str(int(-1)))
256
call my_add_attribute(xf,"spin","-1")
258
call xml_NewElement(xf,"radfunc")
259
call xml_NewElement(xf,"grid")
260
call my_add_attribute(xf,"type",gridtype)
261
call my_add_attribute(xf,"units",gridunits)
262
call my_add_attribute(xf,"scale",gridscale)
263
call my_add_attribute(xf,"step",gridstep)
264
call my_add_attribute(xf,"npts",gridnpoint)
265
call xml_EndElement(xf,"grid")
267
call xml_NewElement(xf,"data")
268
call xml_AddCharacters(xf,viod(ivps,2:nr))
269
call xml_EndElement(xf,"data")
270
call xml_EndElement(xf,"radfunc")
271
call xml_EndElement(xf,"vps")
274
! Up pseudopotentials follows
276
vpsu: do ivps = 1, lmax
277
if (indu(ivps) .eq. 0) cycle
278
call xml_NewElement(xf,"vps")
279
call my_add_attribute(xf,"principal-n",str(no(indu(ivps))))
280
call my_add_attribute(xf,"l",il(ivps))
281
call my_add_attribute(xf,"cutoff",str(rc(ivps)))
282
call my_add_attribute(xf,"occupation",str(zo(indu(ivps))))
283
! JMS/AG 2009/04/02: to be checked which one works better
284
! call my_add_attribute(xf,"spin",str(int(+1)))
285
call my_add_attribute(xf,"spin","+1")
287
call xml_NewElement(xf,"radfunc")
288
call xml_NewElement(xf,"grid")
289
call my_add_attribute(xf,"type",gridtype)
290
call my_add_attribute(xf,"units",gridunits)
291
call my_add_attribute(xf,"scale",gridscale)
292
call my_add_attribute(xf,"step",gridstep)
293
call my_add_attribute(xf,"npts",gridnpoint)
294
call xml_EndElement(xf,"grid")
296
call xml_NewElement(xf,"data")
297
call xml_AddCharacters(xf,viou(ivps,2:nr))
298
call xml_EndElement(xf,"data")
299
call xml_EndElement(xf,"radfunc")
300
call xml_EndElement(xf,"vps")
302
call xml_EndElement(xf,"semilocal")
304
! Dump of the pseudowave functions
305
call xml_NewElement(xf,"pseudowave-functions")
306
call my_add_attribute(xf,"units",slwfnunits)
307
call my_add_attribute(xf,"format",slwfnformat)
308
call my_add_attribute(xf,"n-pseudowave-functions-down",
310
call my_add_attribute(xf,"n-pseudowave-functions-up",
313
! Down pseudowave function follows
315
pswfd: do ivps = 1, lmax
316
if (indd(ivps) .eq. 0) cycle
317
call xml_NewElement(xf,"pswf")
318
call my_add_attribute(xf,"principal-n",str(no(indd(ivps))))
319
call my_add_attribute(xf,"l",il(ivps))
320
call my_add_attribute(xf,"spin","-1")
322
call xml_NewElement(xf,"radfunc")
323
call xml_NewElement(xf,"grid")
324
call my_add_attribute(xf,"type",gridtype)
325
call my_add_attribute(xf,"units",gridunits)
326
call my_add_attribute(xf,"scale",gridscale)
327
call my_add_attribute(xf,"step",gridstep)
328
call my_add_attribute(xf,"npts",gridnpoint)
329
call xml_EndElement(xf,"grid")
331
call xml_NewElement(xf,"data")
332
call xml_AddCharacters(xf,pswfnrd(ivps,2:nr))
333
call xml_EndElement(xf,"data")
334
call xml_EndElement(xf,"radfunc")
335
call xml_EndElement(xf,"pswf")
338
! Up pseudowavefunction follows
340
pswfu: do ivps = 1, lmax
341
if (indu(ivps) .eq. 0) cycle
342
call xml_NewElement(xf,"pswf")
343
call my_add_attribute(xf,"principal-n",str(no(indu(ivps))))
344
call my_add_attribute(xf,"l",il(ivps))
345
call my_add_attribute(xf,"spin","+1")
347
call xml_NewElement(xf,"radfunc")
348
call xml_NewElement(xf,"grid")
349
call my_add_attribute(xf,"type",gridtype)
350
call my_add_attribute(xf,"units",gridunits)
351
call my_add_attribute(xf,"scale",gridscale)
352
call my_add_attribute(xf,"step",gridstep)
353
call my_add_attribute(xf,"npts",gridnpoint)
354
call xml_EndElement(xf,"grid")
356
call xml_NewElement(xf,"data")
357
call xml_AddCharacters(xf,pswfnru(ivps,2:nr))
358
call xml_EndElement(xf,"data")
359
call xml_EndElement(xf,"radfunc")
360
call xml_EndElement(xf,"pswf")
362
call xml_EndElement(xf,"pseudowave-functions")
364
call xml_NewElement(xf,"valence-charge")
365
call xml_NewElement(xf,"radfunc")
366
call xml_NewElement(xf,"grid")
367
call my_add_attribute(xf,"type",gridtype)
368
call my_add_attribute(xf,"units",gridunits)
369
call my_add_attribute(xf,"scale",gridscale)
370
call my_add_attribute(xf,"step",gridstep)
371
call my_add_attribute(xf,"npts",gridnpoint)
372
call xml_EndElement(xf,"grid")
374
call xml_NewElement(xf,"data")
375
call xml_AddCharacters(xf,chval(2:nr))
376
call xml_EndElement(xf,"data")
377
call xml_EndElement(xf,"radfunc")
378
call xml_EndElement(xf,"valence-charge")
380
call xml_NewElement(xf,"pseudocore-charge")
381
call xml_NewElement(xf,"radfunc")
382
call xml_NewElement(xf,"grid")
383
call my_add_attribute(xf,"type",gridtype)
384
call my_add_attribute(xf,"units",gridunits)
385
call my_add_attribute(xf,"scale",gridscale)
386
call my_add_attribute(xf,"step",gridstep)
387
call my_add_attribute(xf,"npts",gridnpoint)
388
call xml_EndElement(xf,"grid")
390
call xml_NewElement(xf,"data")
391
call xml_AddCharacters(xf,cdc(2:nr))
392
call xml_EndElement(xf,"data")
393
call xml_EndElement(xf,"radfunc")
394
call xml_EndElement(xf,"pseudocore-charge")
396
call xml_EndElement(xf,"pseudo")
403
subroutine my_add_attribute(xf,name,value)
404
type(xmlf_t), intent(inout) :: xf
405
character(len=*), intent(in) :: name
406
character(len=*), intent(in) :: value
408
call xml_AddAttribute(xf,name,trim(value))
409
end subroutine my_add_attribute
411
end subroutine pseudoXML