MAISE library interface
General interface function
MAISE functions can be called from external Fortan, C, or C++
codes. The main routine that can be called from external codes is
CALL_MAISE(). This is an interface function that returns the enthalpy,
forces, and stresses calculated with the neural network (NN) or
traditional interatomic potential specified in the model
file. The
input and output variables are explained below, followed by examples
of how the function can be called from basic Fortran, C, and C++
codes.
//======================================================
// Interface for calculating energy, forces, stresses
//------------------------------------------------------
// pointers to be passed from external code
//------------------------------------------------------
// ANN *R neural network model
// PRS *P parser
// PRS *W array of parsers
// LNK *L link
// Cell *C cell
//------------------------------------------------------
// job/cell settings to be defined by user
//------------------------------------------------------
// int CODE external code type (for future use)
// int N number of atoms
// int NM max number of nearest neighbors
// int ND clusters (0) or crystals (3)
// int NP number of cores for parallelization
// int XT fractional (0) or Cartesian (1)
// int *ATMN species types
// double *LAT lattice in 1D array
// double *X positions in 1D array
//------------------------------------------------------
// calculated values
//------------------------------------------------------
// double FRC forces
// double STR stresses
//======================================================
double CALL_MAISE(ANN *R, PRS *P, PRS *W, LNK *L, Cell *C,
int CODE, int N, int NM, int ND, int NP, int XT, int *ATMN, double *LAT, double *X,
double *FRC, double *STR)
{
int i, q;
double H;
C->N = N;
C->NM = NM;
C->ND = ND;
C->NP = NP;
C->XT = XT;
C->R0 = 1.0;
R->JOBT = C->JOBT = 20;
R->NB = C->NB = 1;
R->A = C->A = C->N;
R->NP = C->NP;
C->p /= eV2GPa;
//===== allocate and initialize everything when CALL_MAISE is called the first time =====
if(L->B == 0)
{
P->LM = P->GM = 1;
READ_MODEL(R,P,C);
Build_ANN(R);
if(R->MODT==1)
READ_ANN(R);
else
READ_POT(C,".");
Build_Cell(C,1);
C->NSPC = C->nspc = P->NSPC = R->NSPC;
for(i=0;i<C->nspc;i++)
C->spcz[i] = C->SPCZ[i] = P->SPCZ[i] = R->SPCZ[i];
Build_LNK(L,C->N,C->NM,P->D,3);
for(i=0;i<C->N;i++)
L->MRK[i] = 1;
if(R->MODT==1)
Build_PRS(P,W,0);
}
//===== convert external code cell to MAISE cell =====
CAST_CELL(C, ATMN, LAT, X);
//===== calculate enthalpy, forces, and stresses =====
H = CELL_FRC(R, P, W, C, L);
//===== copy forces and stresses to designated arrays =====
for(i=0;i<C->N;i++)
for(q=0;q<D3;q++)
FRC[3*i+q] = C->F[i][q];
for(q=0;q<6;q++)
STR[q] = C->U[q];
return H;
}
//======================================================
Interface program examples
Examples of basic Fortran, C, and CPP interface programs can be found here.
The interface code will call MAISE to calculate force components and energy for a small nanoparticle using the provided NN model. In order to compile and test the interface program for the desired language (xx = f90, c, or cpp), you should:
download and install MAISE, which will save
libmaise.a
in the installation directorygo to the directory for the desired language (xx = f90, c, or cpp)
set ‘M_LIB’ in the makefile to the directory with installed MAISE
compile the interface code by running
make
go to
./bin
and runxx2maise
.
Note
The initialization of the NN is quite expensive, so it is done only once upon the first CALL_MAISE() call. For this reason, the relevant MAISE structures (ANN, LNK, etc.) are declared in the external code, so that they can be reused in subsequent calls.
LAMMPS interface
MAISE-LAMMPS adds the MAISE-based NN interaction description to the
list of LAMMPS potential types. Once the interface is installed, a
user can run LAMMPS standard simulation jobs (molecular dynamics,
local optimization, etc.) during which MAISE is called to evaluate the
interatomic interactions (energy/enthalpy, forces, and stresses) using
a NN specified in model
.
Installation
If you already have LAMMPS, MAISE-LAMMPS can be installed as follows:
cd lammps/src
git clone https://github.com/maise-guide/maise-lammps.git MAISE; make -C MAISE
If you do not have LAMMPS, please follow these steps:
git clone https://github.com/lammps/lammps.git
cd lammps/src
git clone https://github.com/maise-guide/maise-lammps.git MAISE
make yes-maise
make mpi --jobs
If you wish to remove MAISE from LAMMPS and remove all MAISE files, please run:
cd lammps/src
make no-MAISE
rm -rf MAISE
Setup
A baseline knowledge of LAMMPS is assumed for users of MAISE-LAMMPS.
As usual, the user needs a LAMMPS data file defining the structure and
a LAMMPS input script specifying the simulation. In addition, the user
needs a MAISE model
file defining the specific interaction
description (a NN or a traditional potential).
MAISE model
A MAISE model
can be obtained from the models/
directory in
the MAISE distribution package. Typical model
files are provided
in the examples directories below.
LAMMPS data file
A standard LAMMPS data file should be used. The only requirement is
that atoms appear in the increasing order of their atomic numbers as
in the model
file. Here is an example of a periodic Sn diamond
structure:
# Converted from POSCAR to lammps format
8 atoms
2 atom types
0.000000 6.652395 xlo xhi
0.000000 6.652395 ylo yhi
0.000000 6.652395 zlo zhi
0.000000 0.000000 0.000000 xy xz yz
Atoms
1 2 0.000000 0.000000 0.000000
2 2 3.326198 3.326198 0.000000
3 2 0.000000 3.326198 3.326198
4 2 3.326198 0.000000 3.326198
5 2 1.663099 1.663099 1.663099
6 2 4.989296 4.989296 1.663099
7 2 1.663099 4.989296 4.989296
8 2 4.989296 1.663099 4.989296
Note
The number of atom types is set to 2 because the NN
model
in this example corresponds to the Li-Sn binary
system. The structure is pure Sn because all atom types are
set to 2. A similar convention should be used for systems
with more species.
LAMMPS input file
Detailed information about LAMMPS standard simulation settings can be found on LAMMPS wiki. The table below gives the full list of commands that may have specific settings to define MAISE-LAMMPS runs.
COMMAND |
Short Description |
---|---|
sets units (must be set to metal) |
|
sets boundary conditions (either f f f or p p p) |
|
declares the atom style LAMMPS will use (must be set to atomic) |
|
sets the masses of each species (in the same order as in the model file) |
|
declares the pair style for the run of LAMMPS (must be set to maise) |
|
declares the compute style (must be maisestress) |
|
modifies thermo to use MAISE stresses (see syntax below) |
|
modifies fix to use MAISE stresses (see syntax below) |
units sets the units for the run of LAMMPS. It should always be set to metal for a run of MAISE-LAMMPS. Namely:
units metal
boundary sets either periodic or fixed boundaries. Currently, MAISE treats structures as either nonperiodic or periodic in all three dimensions. So, if a structure has any periodic boundaries then all boundaries must be set to periodic:
boundary p p p
atom_style sets atom style for the run of LAMMPS. It must be atomic for MAISE-LAMMPS:
atom_style atomic
mass sets the masses for all species in the same order they appear in the MAISE model
file. For example:
mass 1 6.941
mass 2 118.71
pair_style declares the pair style for the run of LAMMPS. The
first argument must be maise in order to use MAISE-LAMMPS. The
second argument must be the cutoff distance for each atom in
angstroms, which can be found at the end of the model
file (see
here for more
details). The third argument is the number of cores for
parallelization, which should be a factor of the total number of
atoms. For example:
pair_style maise 7.5 4
compute declares a compute style for a particular value for the run of LAMMPS. In MAISE-LAMMPS, the compute style maisestress is used to pass the stress tensor from MAISE to LAMMPS. The first argument is the compute ID which will be used in the modify commands.
compute press_ID all maisestress thermo_temp
fix_modify modifies a fix to use MAISE stresses. For example:
fix_modify fix_ID energy no virial no press press_ID
Please see the example scripts in src/MAISE/examples/
for basic structure relaxations.
Execution
The standard command is used to run MAISE-LAMMPS jobs. The user should
either make the lmp_mpi
executable globally available or use a
specific path to the executable to run it. An input script
file,
such as in.crystal
or in.cluster
shown below, should be passed
to LAMMPS as follows:
lmp_mpi -in script
Examples
Crystal relaxation
An example of a full crystal structure relaxation at 10 GPa (100 kB)
is given in MAISE/examples/relax-crystal
, with the following
in.crystal
input file:
#MAISE-LAMMPS Example Script 1
units metal
boundary p p p
atom_style atomic
read_data data.crystal
compute press_ID all maisestress thermo_temp
mass 1 6.941
mass 2 118.71
pair_style maise 7.5 4
dump dump_ID all custom 25 exampledump id element x y z fx fy fz
thermo 1
thermo_style custom step ecouple enthalpy press pe
thermo_modify press press_ID
fix fix_ID all box/relax tri 100000.0 vmax .01
fix_modify fix_ID energy no virial no press press_ID
min_style cg
minimize 1.0e-8 1.0e-8 10000 10000
The job should be run in the example directory:
cd lammps/src; cd MAISE/examples/relax-crystal
../../../lmp_mpi -in in.crystal
Note
If a crystal minimization stops before reaching the energy tolerance and has the stopping criteria “linesearch alpha is zero”, tweaking the value vmax in fix box/relax may help. This is not a MAISE-LAMMPS exclusive fix, but it is worthy of note.
Cluster relaxation
An example of a nonperiodic cluster relaxation is given in
MAISE/examples/relax-cluster
, with the following in.cluster
input file:
#MAISE-LAMMPS Example Script 2 (nano)
units metal
boundary f f f
atom_style atomic
read_data data.cluster
mass 1 107.8682
mass 2 196.96657
pair_style maise 7.5 4
dump dump_ID all custom 25 exampledump id element x y z fx fy fz
thermo 1
thermo_style custom step pe
min_style cg
minimize 1.0e-8 1.0e-8 10000 10000
The job should be run in the example directory:
cd lammps/src; cd MAISE/examples/relax-cluster
../../../lmp_mpi -in in.cluster
Note
As usual, fix box/relax is not required for nonperiodic structure relaxations.