TABLE OF CONTENTS


ABINIT/m_ipi [ Modules ]

[ Top ] [ Modules ]

NAME

  m_ipi

FUNCTION

  This module implements the client-side parf of the i-pi protocol.
  The communication between the driver and the client (Abinit) is implemented as follows:
     1)
     2)

  Only structural relaxations with ionmov 28 are implemented.
  Support MD runs requires the proper handling of velocities.
  No parallelism over images.

COPYRIGHT

  Copyright (C) 2020-2024 ABINIT group (MG)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

NOTES

  From the i-pi user guide available at: https://github.com/i-pi/i-pi/blob/master/doc/manual.tex

  The server assumes that 4-byte integers, 8-byte floats and 1-byte characters are used.
  The typical communication flow is as follows:

  1. a header string "STATUS" is sent by the server to the client that has connected to it;

  2. a header string is then returned, giving the status of the client code.
     Recognized messages are:

  "NEEDINIT": if the client code needs any initialising data, it can be sent here.
  The server code will then send a header string "INIT", followed by an integer
  corresponding to the bead index, another integer giving the number of bits
  in the initialization string, and finally the initialization string itself.

  "READY": sent if the client code is ready to calculate the forces. The server
  socket will then send a string "POSDATA", then nine floats for the cell vector
  matrix, then another nine floats for the inverse matrix. The server socket
  will then send one integer giving the number of atoms, then the position data
  as 3 floats for each atom giving the 3 cartesian components of its position.

  "HAVEDATA": is sent if the client has finished computing the potential and
  forces. The server socket then sends a string "GETFORCE", and the client
  socket returns "FORCEREADY". The potential is then returned as a float,
  the number of atoms as an integer, then the force data as 3 floats per atom
  in the same way as the positions, and the virial as 9 floats in the same
  way as the cell vector matrix. Finally, the client may return an arbitrary
  string containing additional data that have been obtained by the electronic
  structure calculation (atomic charges, dipole moment). The client first
  returns an integer specifying the number of characters, and then the string,
  which will be output verbatim if this "extra" information is requested in the
  output section (see 3.2.2).

  3. The server socket waits until the force data for each replica of the system has
  been calculated and returned, then the MD can be propagated for one more time
  step, and new force requests will be dispatched.

SOURCE

62 #if defined HAVE_CONFIG_H
63 #include "config.h"
64 #endif
65 
66 #include "abi_common.h"
67 
68 module m_ipi
69 
70  use defs_basis
71  use m_abicore
72  use m_abimover
73  use m_abihist
74  use m_xmpi
75  use m_errors
76  use m_fsockets
77 
78  use m_fstrings,  only : sjoin, itoa
79  use m_geometry,  only : xcart2xred, det3r, stress_voigt_to_mat
80 
81  implicit none
82 
83  private
84 
85  public :: ipi_setup
86  public :: ipi_shutdown
87  public :: ipi_pred
88  public :: ipi_check_initial_consistency

m_ipi/handle_posdata [ Functions ]

[ Top ] [ m_ipi ] [ Functions ]

NAME

 handle_posdata

FUNCTION

 Receive new geometry from the server after POSDATA header.

INPUTS

OUTPUT

OUTPUT

SOURCE

239 subroutine handle_posdata(out_natom, out_rprimd, out_xred)
240 
241  integer,intent(out) :: out_natom
242  real(dp),intent(out) :: out_rprimd(3,3)
243  real(dp),allocatable,intent(out) :: out_xred(:,:)
244 
245  real(dp) :: gprimd(3,3), mtxbuffer(9)
246  real(dp), allocatable :: combuf(:), xcart(:,:)
247 
248  call readbuffer(socket, mtxbuffer, 9)
249  out_rprimd = transpose(reshape(mtxbuffer, [3, 3]))
250  call readbuffer(socket, mtxbuffer, 9)
251  gprimd = transpose(reshape(mtxbuffer, [3, 3]))
252  call readbuffer(socket, out_natom)
253  ABI_MALLOC(combuf, (3 * out_natom))
254  call readbuffer(socket, combuf, 3 * out_natom)
255  ABI_MALLOC(xcart, (3, out_natom))
256  xcart = reshape(combuf, [3, out_natom])
257  ABI_MALLOC(out_xred, (3, out_natom))
258  call xcart2xred(out_natom, out_rprimd, xcart, out_xred)
259  ABI_FREE(combuf)
260  ABI_FREE(xcart)
261 
262 end subroutine handle_posdata

m_ipi/ipi_check_initial_consistency [ Functions ]

[ Top ] [ m_ipi ] [ Functions ]

NAME

 ipi_check_initial_consistency

FUNCTION

 Compare initial configuration read from input file with the one trasmitted by the server.

INPUTS

OUTPUT

OUTPUT

SOURCE

177 subroutine ipi_check_initial_consistency(in_natom, in_rprimd, in_xred, ierr)
178 
179  integer,intent(in) :: in_natom
180  real(dp),intent(in) :: in_rprimd(3,3)
181  real(dp),intent(in) :: in_xred(3,in_natom)
182  integer,intent(out) :: ierr
183 
184  integer :: ii
185 
186  call wrtout(std_out, "ipi mode: Checking whether initial geometry from server agrees with input file")
187 
188  ierr = 0
189  if (in_natom /= origin_natom) then
190    ABI_WARNING(sjoin("in_natom:", sjoin(itoa(in_natom), " != origin_natom", itoa(origin_natom))))
191    ierr = ierr + 1
192  end if
193 
194  if (any(abs(in_rprimd - origin_rprimd) > tol6)) then
195    ABI_WARNING("Mismatch between input file and data from socket: in_rprimd and origin_rprimd do not agree within 1e-6")
196    write(std_out, "(a)")" in_rprind(:,ii), origin_rprimd(:,ii)"
197    do ii=1,3
198      write(std_out, *)in_rprimd(:,ii), origin_rprimd(:,ii)
199    end do
200    ierr = ierr + 1
201  end if
202 
203  if (.not. allocated(origin_xred)) then
204    ierr = ierr + 1
205    ABI_WARNING("origin_xred is not allocated!")
206  end if
207 
208  if (in_natom == origin_natom .and. allocated(origin_xred)) then
209    if (any(abs(in_xred - origin_xred) > tol6)) then
210      ABI_WARNING("Mismatch between input file and data from socket: in_xred and origin_xred do not agree withing 1e-6")
211      ierr = ierr + 1
212      write(std_out, "(a)")" in_xred(:,ii), origin_xred(:,ii)"
213      do ii=1,3
214        write(std_out, *)in_xred(:,ii), origin_xred(:,ii)
215      end do
216    end if
217  end if
218 
219  !ABI_CHECK(ierr == 0, "Initial structure sent by i-pi server does not match the one read from file! See messages above")
220 
221 end subroutine ipi_check_initial_consistency

m_ipi/ipi_pred [ Functions ]

[ Top ] [ m_ipi ] [ Functions ]

NAME

 ipi_pred

FUNCTION

  "predict" new structure using the data sent by the server.

INPUTS

 ab_mover <type(abimover)> : Datatype with all the information needed by the predictor
 itime  : Index of the present iteration
 ntime  : Maximal number of iterations
 zDEBUG : if true print some debugging information

OUTPUT

SIDE EFFECTS

 hist <type(abihist)>: History of positions,forces,acell, rprimd, stresses

OUTPUT

SOURCE

313 subroutine ipi_pred(ab_mover, hist, itime, ntime, zDEBUG, iexit, comm_cell)
314 
315 !Arguments ------------------------------------
316 !scalars
317  integer,intent(in) :: itime, ntime, iexit, comm_cell
318  logical,intent(in) :: zDEBUG
319  type(abimover),intent(in) :: ab_mover
320  type(abihist),intent(inout) :: hist
321 
322 !local variables
323  integer,parameter :: master = 0
324  integer :: my_rank, natom, new_natom, ierr !, ihist_prev
325  real(dp) :: etotal, ucvol
326  character(len=HDRLEN) :: header
327  character(len=500) :: msg
328 !arrays
329  real(dp) :: acell(3), rprimd(3,3), xred(3,ab_mover%natom), strten(6), sigma(3,3)
330  real(dp) :: new_rprimd(3, 3)
331  real(dp),allocatable :: new_xred(:,:)
332 
333 ! *************************************************************************
334 
335  write(msg, "(3(a,i0))")" ipi_pred: itime: ", itime, ", ntime: ", ntime, ", iexit:", iexit
336  call wrtout(std_out, msg)
337 
338  natom = ab_mover%natom
339  ! Obtain the present values from the history
340  call hist2var(acell, hist, natom, rprimd, xred, zDEBUG)
341  etotal = hist%etot(hist%ihist)
342  strten = hist%strten(:,hist%ihist)
343  ucvol = det3r(rprimd)
344 
345  if (my_rank == master) then
346    call readbuffer(socket, header, HDRLEN)
347    ABI_CHECK(trim(header) == "STATUS", sjoin("Expecting STATUS header, got:", trim(header)))
348    call writebuffer( socket, "HAVEDATA    ", HDRLEN )
349    call readbuffer(socket, header, HDRLEN)
350    ABI_CHECK(trim(header) == "GETFORCE", sjoin("Expecting GETFORCE header, got:", trim(header)))
351 
352    ! Communicate energy info back to i-pi
353    call wrtout(std_out, " i-pi mode: Returning etotal, forces, and stress tensor to server...")
354    call writebuffer(socket, "FORCEREADY  ", HDRLEN)
355    call writebuffer(socket, etotal)
356    call writebuffer(socket, natom)
357    call writebuffer_dv(socket, hist%fcart(:,:,hist%ihist), 3 * natom)
358    call stress_voigt_to_mat(strten, sigma)
359    ! Well it's a symmetric tensor.
360    ! TODO: Check volume
361    sigma = -ucvol * transpose(sigma)
362    call writebuffer_dv(socket, sigma, 9)
363 
364    ! Note: i-pi can also receive an arbitrary string, that will be printed
365    ! out to the "extra" trajectory file. This is useful if you want to
366    ! return additional information, e.g. atomic charges, wannier centres,
367    ! etc. one must return the number of characters, then the string. Here
368    ! we just send back zero characters.
369    !
370    call writebuffer(socket, 0)
371 
372    ! ===============================================================
373    ! Here master waits for new lattice and positions from the server
374    ! ===============================================================
375    call readbuffer(socket, header, HDRLEN)
376    ABI_CHECK(trim(header) == "STATUS", sjoin("Expecting STATUS header, got:", trim(header)))
377    call writebuffer(socket, "READY       ", HDRLEN)
378    call readbuffer(socket, header, HDRLEN)
379    ABI_CHECK(trim(header) == "POSDATA", sjoin("Expecting POSDATA header, got:", trim(header)))
380    call handle_posdata(new_natom, new_rprimd, new_xred)
381 
382    ! Some basic consistency checks. Obviously there are lot of things that can go wrong here
383    ! if the driver changes the space group or the input file is not consistent with
384    ! the logic implemented by the server.
385    ABI_CHECK(new_natom == ab_mover%natom, "ipi server shall not change the number of atoms!")
386    if (ab_mover%optcell == 0 .and. any(abs(new_rprimd - origin_rprimd) > tol6)) then
387      ABI_ERROR("Mismatch between origin_rprimd and data from socket: origin_rprimd and new_rprimd do not agree within 1e-6")
388    end if
389  end if
390 
391  ! ====================
392  ! All procs block here
393  ! ====================
394  call xmpi_barrier(comm_cell)
395  call xmpi_bcast(new_rprimd, master, comm_cell, ierr)
396  call xmpi_bcast(new_xred, master, comm_cell, ierr)
397 
398  ! Update the history with the prediction
399  ! Increase indexes
400  hist%ihist = abihist_findIndex(hist, +1)
401  !call mkradim(acell,rprim,rprimd_symm)
402 
403  ! Fill the history with the variables xred, acell, rprimd, vel
404  acell = one
405  call var2hist(acell, hist, ab_mover%natom, new_rprimd, new_xred, zDEBUG)
406  !ihist_prev = abihist_findIndex(hist, -1)
407  ! FIXME: I don't know how to handle vel if MD run!
408  !hist%vel(:,:,hist%ihist) = hist%vel(:,:,ihist_prev)
409 
410 end subroutine ipi_pred

m_ipi/ipi_setup [ Functions ]

[ Top ] [ m_ipi ] [ Functions ]

NAME

 ipi_setup

FUNCTION

  Initialize the socket from a string that is usually passed via the command line interface.
  Get the initial configuration from the server and save it in module global variables
  for subsequent consistency check.
  These values, indeed, must agree with those read from the input file.
  See also ipi_check_initial_consistency.

INPUTS

OUTPUT

OUTPUT

SOURCE

123 subroutine ipi_setup(string, comm)
124 
125 !Arguments ------------------------------------
126 !scalars
127  character(len=*),intent(in) :: string
128  integer,intent(in) :: comm
129 
130 !local variables
131  integer,parameter :: master = 0
132  integer :: my_rank, ierr
133  character(len=HDRLEN) :: header
134 
135 ! *************************************************************************
136 
137  my_rank = xmpi_comm_rank(comm)
138 
139  call wrtout(std_out, "Initializing i-pi protocol...")
140  if (my_rank == master) then
141    call socket_from_string(string, socket)
142    call readbuffer(socket, header, HDRLEN)
143    ABI_CHECK(trim(header) == "STATUS", sjoin("Expecting STATUS header, got:", trim(header)))
144    call writebuffer(socket, "READY       ", HDRLEN)
145    call readbuffer(socket, header, HDRLEN)
146    ABI_CHECK(trim(header) == "POSDATA", sjoin("Expecting POSDATA header, got:", trim(header)))
147    call handle_posdata(origin_natom, origin_rprimd, origin_xred)
148  end if
149 
150  ! broadcast data to other procs.
151  call xmpi_bcast(origin_natom, master, comm, ierr)
152  call xmpi_bcast(origin_rprimd, master, comm, ierr)
153  if (my_rank /= master) then
154    ABI_MALLOC(origin_xred, (3, origin_natom))
155  end if
156  call xmpi_bcast(origin_xred, master, comm, ierr)
157  call wrtout(std_out, "ipi_setup completed")
158 
159 end subroutine ipi_setup

m_ipi/ipi_shutdown [ Functions ]

[ Top ] [ m_ipi ] [ Functions ]

NAME

 ipi_shutdown

FUNCTION

INPUTS

OUTPUT

OUTPUT

SOURCE

279 subroutine ipi_shutdown()
280 
281 ! *************************************************************************
282 
283  origin_natom = 0
284  origin_rprimd = zero
285  ABI_SFREE(origin_xred)
286  !call close_socket()
287 
288 end subroutine ipi_shutdown