update examples/COUPLE/simple sources to work with the current version of LAMMPS

This commit is contained in:
Axel Kohlmeyer
2020-09-24 01:30:34 -04:00
parent b350dce573
commit ac78f8f1e1
5 changed files with 40 additions and 40 deletions

View File

@ -2,3 +2,7 @@
*.d *.d
*.a *.a
*~ *~
liblammps.mod
simpleC
simpleCC
simpleF

View File

@ -8,16 +8,15 @@ code to perform a coupled calculation.
simple.cpp is the C++ driver simple.cpp is the C++ driver
simple.c is the C driver simple.c is the C driver
simple.f90 is the Fortran driver simple.f90 is the Fortran driver
libfwrapper.c is the Fortran-to-C wrapper
The 3 codes do the same thing, so you can compare them to see how to The 3 codes do the same thing, so you can compare them to see how to
drive LAMMPS from each language. See lammps/python/example/simple.py drive LAMMPS from each language. See python/example/simple.py
to do something similar from Python. The Fortran driver requires an to do something similar from Python. The Fortran driver requires an
additional wrapper library that interfaces the C interface of the additional wrapper library that interfaces the C interface of the
LAMMPS library to Fortran and also translates the MPI communicator LAMMPS library to Fortran and also translates the MPI communicator
from Fortran to C. from Fortran to C.
First build LAMMPS as a library (see examples/COUPLE/README), e.g. First build LAMMPS as a library (see examples/COUPLE/README), e.g.
make mode=shlib mpi make mode=shlib mpi
@ -26,27 +25,23 @@ these, which include paths to the LAMMPS library interface, and
linking with FFTW (only needed if you built LAMMPS as a library with linking with FFTW (only needed if you built LAMMPS as a library with
its PPPM solver). its PPPM solver).
This builds the C++ driver with the LAMMPS library using the mpiCC This builds the C++ driver with the LAMMPS library using the mpicxx
(C++) compiler: (C++) compiler:
mpiCC -I/home/sjplimp/lammps/src -c simple.cpp mpicxx -I/home/sjplimp/lammps/src -c simple.cpp
mpiCC -L/home/sjplimp/lammps/src simple.o -llammps -lfftw -o simpleCC mpicxx -L/home/sjplimp/lammps/src simple.o -llammps -o simpleCC
This builds the C driver with the LAMMPS library using the mpicc (C) This builds the C driver with the LAMMPS library using the mpicc (C)
compiler: compiler:
mpicc -I/home/sjplimp/lammps/src -c simple.c mpicc -I/home/sjplimp/lammps/src -c simple.c
mpicc -L/home/sjplimp/lammps/src simple.o -llammps -lfftw -o simpleC mpicc -L/home/sjplimp/lammps/src simple.o -llammps -o simpleC
This builds the Fortran wrapper and driver with the LAMMPS library This builds the Fortran module and driver with the LAMMPS library
using the mpicc (C) and mpifort (Fortran) compilers, using the wrapper using the mpifort (Fortran) compilers, using the Fortran module from
in the fortran directory: the fortran directory:
cp ../fortran/libfwrapper.c . mpifort -L/home/sjplimp/lammps/src ../../../fortran/lammps.f90 simple.f90 -llammps -o simpleF
mpicc -I/home/sjplimp/lammps/src -c libfwrapper.c
mpifort -c simple.f90
mpifort -L/home/sjplimp/lammps/src simple.o libfwrapper.o \
-llammps -lfftw -o simpleF
You then run simpleCC, simpleC, or simpleF on a parallel machine You then run simpleCC, simpleC, or simpleF on a parallel machine
on some number of processors Q with 2 arguments: on some number of processors Q with 2 arguments:

View File

@ -23,7 +23,7 @@
#include <stdlib.h> #include <stdlib.h>
#include <string.h> #include <string.h>
#include <mpi.h> #include <mpi.h>
#include "lammps/library.h" /* this is a LAMMPS include file */ #include "library.h" /* this is a LAMMPS include file */
int main(int narg, char **arg) int main(int narg, char **arg)
{ {
@ -124,10 +124,10 @@ int main(int narg, char **arg)
/* use commands_string() and commands_list() to invoke more commands */ /* use commands_string() and commands_list() to invoke more commands */
char *strtwo = "run 10\nrun 20"; const char *strtwo = "run 10\nrun 20";
if (lammps == 1) lammps_commands_string(lmp,strtwo); if (lammps == 1) lammps_commands_string(lmp,strtwo);
char *cmds[2]; const char *cmds[2];
cmds[0] = "run 10"; cmds[0] = "run 10";
cmds[1] = "run 20"; cmds[1] = "run 20";
if (lammps == 1) lammps_commands_list(lmp,2,cmds); if (lammps == 1) lammps_commands_list(lmp,2,cmds);

View File

@ -25,10 +25,10 @@
#include <mpi.h> #include <mpi.h>
// these are LAMMPS include files // these are LAMMPS include files
#include <lammps/lammps.h> #include "lammps.h"
#include <lammps/input.h> #include "input.h"
#include <lammps/atom.h> #include "atom.h"
#include <lammps/library.h> #include "library.h"
using namespace LAMMPS_NS; using namespace LAMMPS_NS;
@ -135,10 +135,10 @@ int main(int narg, char **arg)
// use commands_string() and commands_list() to invoke more commands // use commands_string() and commands_list() to invoke more commands
char *strtwo = (char *) "run 10\nrun 20"; const char *strtwo = (char *) "run 10\nrun 20";
if (lammps == 1) lammps_commands_string(lmp,strtwo); if (lammps == 1) lammps_commands_string(lmp,strtwo);
char *cmds[2]; const char *cmds[2];
cmds[0] = (char *) "run 10"; cmds[0] = (char *) "run 10";
cmds[1] = (char *) "run 20"; cmds[1] = (char *) "run 20";
if (lammps == 1) lammps_commands_list(lmp,2,cmds); if (lammps == 1) lammps_commands_list(lmp,2,cmds);

View File

@ -18,13 +18,14 @@
! See README for compilation instructions ! See README for compilation instructions
PROGRAM f_driver PROGRAM f_driver
USE mpi
USE liblammps
IMPLICIT NONE IMPLICIT NONE
INCLUDE 'mpif.h'
INTEGER, PARAMETER :: fp=20 INTEGER, PARAMETER :: fp=20
INTEGER :: n, narg, ierr, me, nprocs, natoms INTEGER :: n, narg, ierr, me, nprocs, natoms
INTEGER :: lammps, nprocs_lammps, comm_lammps INTEGER :: color, nprocs_lammps, comm_lammps
INTEGER (kind=8) :: ptr TYPE(LAMMPS) :: lmp
REAL (kind=8), ALLOCATABLE :: x(:) REAL (kind=8), ALLOCATABLE :: x(:)
REAL (kind=8), PARAMETER :: epsilon=0.1 REAL (kind=8), PARAMETER :: epsilon=0.1
@ -58,14 +59,14 @@ PROGRAM f_driver
END IF END IF
END IF END IF
lammps = 0 color = 0
IF (me < nprocs_lammps) THEN IF (me < nprocs_lammps) THEN
lammps = 1 color = 1
ELSE ELSE
lammps = MPI_UNDEFINED color = MPI_UNDEFINED
END IF END IF
CALL mpi_comm_split(MPI_COMM_WORLD,lammps,0,comm_lammps,ierr) CALL mpi_comm_split(MPI_COMM_WORLD,color,0,comm_lammps,ierr)
! open LAMMPS input script on rank zero ! open LAMMPS input script on rank zero
@ -81,7 +82,7 @@ PROGRAM f_driver
! (could just send it to proc 0 of comm_lammps and let it Bcast) ! (could just send it to proc 0 of comm_lammps and let it Bcast)
! all LAMMPS procs call lammps_command() on the line */ ! all LAMMPS procs call lammps_command() on the line */
IF (lammps == 1) CALL lammps_open(comm_lammps,ptr) IF (color == 1) lmp=lammps(comm=comm_lammps)
n = 0 n = 0
DO DO
@ -99,7 +100,7 @@ PROGRAM f_driver
CALL mpi_bcast(n,1,MPI_INTEGER,0,MPI_COMM_WORLD,ierr) CALL mpi_bcast(n,1,MPI_INTEGER,0,MPI_COMM_WORLD,ierr)
IF (n == 0) EXIT IF (n == 0) EXIT
CALL mpi_bcast(line,n,MPI_CHARACTER,0,MPI_COMM_WORLD,ierr) CALL mpi_bcast(line,n,MPI_CHARACTER,0,MPI_COMM_WORLD,ierr)
IF (lammps == 1) CALL lammps_command(ptr,line,n) IF (color == 1) CALL lmp%command(line(1:n))
END DO END DO
CLOSE(UNIT=fp) CLOSE(UNIT=fp)
@ -109,27 +110,27 @@ PROGRAM f_driver
! put coords back into LAMMPS ! put coords back into LAMMPS
! run a single step with changed coords */ ! run a single step with changed coords */
IF (lammps == 1) THEN IF (color == 1) THEN
CALL lammps_command(ptr,'run 10',6) CALL lmp%command('run 10')
CALL lammps_get_natoms(ptr,natoms) natoms = NINT(lmp%get_natoms())
ALLOCATE(x(3*natoms)) ALLOCATE(x(3*natoms))
! these calls are commented out, b/c libfwrapper.c ! these calls are commented out, b/c libfwrapper.c
! needs to be updated to use gather_atoms and scatter_atoms ! needs to be updated to use gather_atoms and scatter_atoms
!CALL lammps_gather_atoms(ptr,'x',1,3,x); !CALL lammps_gather_atoms(ptr,'x',1,3,x)
!x(1) = x(1) + epsilon !x(1) = x(1) + epsilon
!CALL lammps_scatter_atoms(ptr,'x',1,3,x); !CALL lammps_scatter_atoms(ptr,'x',1,3,x)
DEALLOCATE(x) DEALLOCATE(x)
CALL lammps_command(ptr,'run 1',5); CALL lmp%command('run 1')
END IF END IF
! free LAMMPS object ! free LAMMPS object
IF (lammps == 1) CALL lammps_close(ptr); IF (color == 1) CALL lmp%close()
! close down MPI ! close down MPI