TABLE OF CONTENTS


ABINIT/out_geometry_xml [ Functions ]

[ Top ] [ Functions ]

NAME

 out_geometry_xml

FUNCTION

 Output in the XML file, the box size and the atomic position.
 (see extras/post_processing/abinitRun.dtd)

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  dtset <type(dataset_type)>=all input variables for this dataset
  level=use for indentation of the XML, two spaces are added by level.
  natom=number of atoms.
  rprimd(3,3)=dimensional primitive translations in real space (bohr)
  xred(3,natom)=reduced dimensionless atomic coordinates

NOTES

PARENTS

      scfcv

CHILDREN

      xred2xcart

SOURCE

324 #if defined HAVE_CONFIG_H
325 #include "config.h"
326 #endif
327 
328 subroutine out_geometry_XML(dtset, level, natom, rprimd, xred)
329 
330  use defs_basis
331  use defs_abitypes
332  use m_profiling_abi
333 
334 !This section has been created automatically by the script Abilint (TD).
335 !Do not modify the following lines by hand.
336 #undef ABI_FUNC
337 #define ABI_FUNC 'out_geometry_XML'
338  use interfaces_41_geometry
339 !End of the abilint section
340 
341  implicit none
342 
343 !Arguments ------------------------------------
344 !scalars
345  integer,intent(in) :: level,natom
346  type(dataset_type),intent(in) :: dtset
347 !arrays
348  real(dp),intent(in) :: rprimd(3,3)
349  real(dp),intent(inout) :: xred(3,natom)
350 
351 !Local variables -------------------------
352   character(len = 1), parameter :: rprimd_names(3) = (/ "x", "y", "z" /)
353 !scalars
354  integer :: i,j
355  character(len=128) :: spacer,value
356 !arrays
357  real(dp),allocatable :: xcart(:,:)
358 
359 ! *************************************************************************
360 
361 !Compute the spacer to put before each markup
362  write(spacer, "(A,I0)") "A", 2 * level
363  write(ab_xml_out, "("//trim(spacer)//",A)") " ", '<geometry>'
364 !Compute the cartesian coordinates of atoms
365  ABI_ALLOCATE(xcart,(3, natom))
366  call xred2xcart(natom, rprimd, xcart, xred)
367 !Ouput the rprimd matrix
368  write(ab_xml_out, "("//trim(spacer)//",A)", advance = "NO") " ", '  <rprimd'
369  do i = 1, 3, 1
370    do j = 1, 3, 1
371      write(value, "(es20.8)") rprimd(i, j)
372      write(ab_xml_out, "(A,A,I0,A,A,A)", advance = "NO") ' ', rprimd_names(i), j, '="', trim(value) ,'"'
373    end do
374  end do
375  write(ab_xml_out, "(A)") ' />'
376 !Output the atom position
377  do i = 1, natom, 1
378    write(ab_xml_out, "("//trim(spacer)//",A)", advance = "NO") " ", '  <position'
379    write(ab_xml_out, "(A,I0,A,I0,A)", advance = "NO") ' atom="a_', dtset%jdtset ,'_', i ,'"'
380    do j = 1, 3, 1
381      write(value, "(es20.8)") xcart(j, i)
382      write(ab_xml_out, "(A,A,A,A,A)", advance = "NO") ' ', rprimd_names(j), '="', trim(value) ,'"'
383    end do
384    write(ab_xml_out, "(A)") ' />'
385  end do
386  ABI_DEALLOCATE(xcart)
387  write(ab_xml_out, "("//trim(spacer)//",A)") " ", '</geometry>'
388 
389 end subroutine out_geometry_XML

ABINIT/out_resultsgs_xml [ Functions ]

[ Top ] [ Functions ]

NAME

 out_resultsgs_xml

FUNCTION

 Output in the XML file, the decomposition of the energy and
 the forces after a scfcv loop.
 (see extras/post_processing/abinitRun.dtd)
 Warning : this method is not thread safe and should be called
 only by one thread.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  dtset <type(dataset_type)>=all input variables for this dataset
  level=use for indentation of the XML, two spaces are added by level.
  results_gs <type(results_gs_type)>=results (energy and its components,
   forces and its components, the stress tensor) of a ground-state computation.
  usepaw= 0 for non paw calculation; =1 for paw calculation

NOTES

PARENTS

      scfcv

CHILDREN

      xred2xcart

SOURCE

187 subroutine out_resultsgs_XML(dtset, level, results_gs, usepaw)
188 
189  use defs_basis
190  use defs_abitypes
191  use m_profiling_abi
192 
193  use m_results_gs , only : results_gs_type
194 
195 !This section has been created automatically by the script Abilint (TD).
196 !Do not modify the following lines by hand.
197 #undef ABI_FUNC
198 #define ABI_FUNC 'out_resultsgs_XML'
199 !End of the abilint section
200 
201  implicit none
202 
203 !Arguments ------------------------------------
204 !scalars
205  integer,intent(in) :: level,usepaw
206  type(dataset_type),intent(in) :: dtset
207  type(results_gs_type),intent(inout) :: results_gs
208 
209 !Local variables -------------------------
210   character(len = 1), parameter :: axes_names(3) = (/ "x", "y", "z" /)
211 !scalars
212  integer :: i,j
213  character(len=128) :: spacer,value
214 
215 ! *************************************************************************
216 
217 !Compute the spacer to put before each markup
218  write(spacer, "(A,I0)") "A", 2 * level
219 
220 !Begin with the energy part
221  write(ab_xml_out, "("//trim(spacer)//",A)", advance = "NO") " ", '<energy'
222  if (dtset%iscf < 10) then
223    write(ab_xml_out, "(A)", advance = "NO") ' type="direct"'
224    write(value, "(es20.8)") results_gs%energies%e_kinetic
225    write(ab_xml_out, "(A,A,A)", advance = "NO") ' kinetic="', trim(value) ,'"'
226    write(value, "(es20.8)") results_gs%energies%e_localpsp
227    write(ab_xml_out, "(A,A,A)", advance = "NO") ' local="', trim(value) ,'"'
228    write(value, "(es20.8)") results_gs%energies%e_nonlocalpsp
229    write(ab_xml_out, "(A,A,A)", advance = "NO") ' non-local="', trim(value) ,'"'
230    if (usepaw == 1) then
231      write(value, "(es20.8)") results_gs%energies%e_paw
232      write(ab_xml_out, "(A,A,A)", advance = "NO") ' paw="', trim(value) ,'"'
233    end if
234  else
235    write(ab_xml_out, "(A)", advance = "NO") ' type="double-counting"'
236    write(value, "(es20.8)") results_gs%energies%e_eigenvalues
237    write(ab_xml_out, "(A,A,A)", advance = "NO") ' eigen-values="', trim(value) ,'"'
238    write(value, "(es20.8)") results_gs%energies%e_xcdc
239    write(ab_xml_out, "(A,A,A)", advance = "NO") ' xcdc="', trim(value) ,'"'
240    if (usepaw == 1) then
241      write(value, "(es20.8)") results_gs%energies%e_pawdc
242      write(ab_xml_out, "(A,A,A)", advance = "NO") ' pawdc="', trim(value) ,'"'
243    end if
244  end if
245  if (dtset%berryopt == 4 .or. dtset%berryopt == 6 .or. dtset%berryopt == 7 .or. &
246 & dtset%berryopt ==14 .or. dtset%berryopt ==16 .or. dtset%berryopt ==17) then
247    write(value, "(es20.8)") results_gs%energies%e_elecfield
248    write(ab_xml_out, "(A,A,A)", advance = "NO") ' electric-field="', trim(value) ,'"'
249  end if
250  if(dtset%occopt >= 3 .and. dtset%occopt <= 8) then
251    write(value, "(es20.8)") results_gs%energies%e_entropy
252    write(ab_xml_out, "(A,A,A)", advance = "NO") ' entropy="', trim(value) ,'"'
253  end if
254  write(value, "(es20.8)") results_gs%energies%e_ewald
255  if (dtset%icoulomb == 0) then
256    write(ab_xml_out, "(A,A,A)", advance = "NO") ' ewald="', trim(value) ,'"'
257  else
258    write(ab_xml_out, "(A,A,A)", advance = "NO") ' ion-ion="', trim(value) ,'"'
259  end if
260  write(value, "(es20.8)") results_gs%energies%e_chempot
261  write(ab_xml_out, "(A,A,A)", advance = "NO") ' chempot="', trim(value) ,'"'
262  write(value, "(es20.8)") results_gs%energies%e_hartree
263  write(ab_xml_out, "(A,A,A)", advance = "NO") ' hartree="', trim(value) ,'"'
264  write(value, "(es20.8)") results_gs%energies%e_corepsp
265  write(ab_xml_out, "(A,A,A)", advance = "NO") ' core="', trim(value) ,'"'
266  write(value, "(es20.8)") results_gs%energies%e_xc
267  write(ab_xml_out, "(A,A,A)", advance = "NO") ' xc="', trim(value) ,'"'
268  write(value, "(es20.8)") results_gs%etotal
269  write(ab_xml_out, "(A,A,A)", advance = "NO") ' total="', trim(value) ,'"'
270  write(ab_xml_out, "(A)") ' />'
271 
272 
273 !finish with the forces part
274  if (dtset%optforces == 1) then
275    write(ab_xml_out, "("//trim(spacer)//",A)") " ", '<forces>'
276    do i = 1, dtset%natom, 1
277      write(ab_xml_out, "("//trim(spacer)//",A)", advance = "NO") " ", '  <force'
278      write(ab_xml_out, "(A,I0,A,I0,A)", advance = "NO") ' atom="a_', dtset%jdtset ,'_', i ,'"'
279      do j = 1, 3, 1
280        write(value, "(es20.8)") results_gs%fcart(j, i)
281        write(ab_xml_out, "(A,A,A,A,A)", advance = "NO") ' ', axes_names(j), '="', trim(value) ,'"'
282      end do
283      write(ab_xml_out, "(A)") ' />'
284    end do
285    write(ab_xml_out, "("//trim(spacer)//",A)") " ", '</forces>'
286  end if
287 
288 end subroutine out_resultsgs_XML

ABINIT/outxml_finalise [ Functions ]

[ Top ] [ Functions ]

NAME

 outxml_finalise

FUNCTION

 Close the XML log file, and write the timing information.
 (see extras/post_processing/abinitRun.dtd)
 Warning : this method is not thread safe and should be called
 only by one thread.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  tsec=the cpu time and the wall time in seconds.
  values=the date values returned by date_and_time() intrinsic Fortran routine.

NOTES

PARENTS

      abinit

CHILDREN

      xred2xcart

SOURCE

111 subroutine outxml_finalise(tsec, values)
112 
113  use defs_basis
114  use m_profiling_abi
115 
116 !This section has been created automatically by the script Abilint (TD).
117 !Do not modify the following lines by hand.
118 #undef ABI_FUNC
119 #define ABI_FUNC 'outxml_finalise'
120 !End of the abilint section
121 
122   implicit none
123 
124 !Arguments -------------------------------
125 
126   integer, intent(in) :: values(8)
127   real(dp), intent(in) :: tsec(2)
128 
129 !Local variables -------------------------
130 
131   character(len=500) :: message
132 
133 ! *********************************************************************
134 
135  write(ab_xml_out, "(A)") '  <!-- timeInfo markup : cpu and wall attributes are given in seconds. -->'
136  write(ab_xml_out, "(A,I0,A,I0,A,I0,A)", advance = "NO") &
137 & '  <timeInfo day="', values(3) , &
138 & '" month="', values(2) ,'" year="', values(1) ,'" '
139  write(message, *) tsec(1)
140  write(ab_xml_out, "(A,A,A)", advance = "NO") 'cpu="', trim(message) ,'"'
141  write(message, *) tsec(2)
142  write(ab_xml_out, "(A,A,A)") ' wall="', trim(message) ,'" />'
143  write(ab_xml_out, "(A)") "</abinitRun>"
144 !ABINIT has been compiled with XML output support, then we close the
145 !channel of the XML output file.
146  close(unit = ab_xml_out)
147 
148 end subroutine outxml_finalise

ABINIT/outxml_open [ Functions ]

[ Top ] [ Functions ]

NAME

 outxml_open

FUNCTION

 Open the XML log file, and write the header inside.
 (see extras/post_processing/abinitRun.dtd)
 Warning : this method is not thread safe and should be called
 only by one thread.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  filename=the name of the file to write to.

NOTES

PARENTS

      abinit

CHILDREN

      xred2xcart

SOURCE

33 #if defined HAVE_CONFIG_H
34 #include "config.h"
35 #endif
36 
37 #include "abi_common.h"
38 
39 subroutine outxml_open(filename)
40 
41  use defs_basis
42  use m_profiling_abi
43  use m_errors
44 
45  use m_io_tools,   only: open_file
46 
47 !This section has been created automatically by the script Abilint (TD).
48 !Do not modify the following lines by hand.
49 #undef ABI_FUNC
50 #define ABI_FUNC 'outxml_open'
51 !End of the abilint section
52 
53   implicit none
54 
55 !Arguments -------------------------------
56 
57   character(len = *), intent(in) :: filename
58 
59 !Local variables -------------------------
60  character(len=500) :: msg
61 
62 ! *********************************************************************
63 
64 !ABINIT has been compiled with XML output support, then we open the
65 !channel of the XML output file.
66  if (open_file(trim(filename)//"_LOG.xml", msg, unit=ab_xml_out, form="formatted", action="write") /= 0) then
67    MSG_ERROR(msg)
68  end if
69 
70  write(ab_xml_out, "(A)") '<?xml version="1.0" encoding="utf-8"?>'
71 !DTD declaration : root element is "abinitRun", and DTD file is not public
72 !but given in the distribution by the file abinitRun.dtd.
73  write(ab_xml_out, "(A)") '<!DOCTYPE abinitRun SYSTEM "extras/post_processing/abinitRun.dtd">'
74 !Creating the root markup.
75  write(ab_xml_out, "(A)") '<abinitRun>'
76 
77 end subroutine outxml_open