Merge pull request #3727 from taylor-a-barnes/taylor

Small enhancements to the mdi/qm fix
This commit is contained in:
Axel Kohlmeyer
2023-04-19 14:22:11 -04:00
committed by GitHub
9 changed files with 438 additions and 300 deletions

View File

@ -294,6 +294,11 @@ LAMMPS can also be used as an MDI driver in other unit choices it
supports, e.g. *lj*, but then no unit conversion to MDI units is
performed.
If this fix is used in conjuction with a QM code that does not support
periodic boundary conditions (more specifically, a QM code that does
not support the ``>CELL`` MDI command), the LAMMPS system must be
fully non-periodic. I.e. no dimension of the system can be periodic.
Related commands
""""""""""""""""

View File

@ -261,6 +261,11 @@ used to specify *real* or *metal* units. This will ensure the correct
unit conversions between LAMMPS and MDI units. The other code will
also perform similar unit conversions into its preferred units.
If this fix is used in conjuction with a QM code that does not support
periodic boundary conditions (more specifically, a QM code that does
not support the ``>CELL`` MDI command), the LAMMPS system must be
fully non-periodic. I.e. no dimension of the system can be periodic.
Related commands
""""""""""""""""

View File

@ -6,7 +6,7 @@ SPPARKS applications in 2 sister directories.
The library dir has a Makefile (which you may need to edit for your
box). If you type
g++ -f Makefile.mpi
make -f Makefile.mpi
you should create libcouple.a, which the other coupled applications
link to.

View File

@ -99,18 +99,14 @@ FixMDIQM::FixMDIQM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
} else if (strcmp(arg[iarg], "elements") == 0) {
const char *symbols[] = {
"H" , "He", "Li", "Be", "B" , "C" , "N" , "O" , "F" , "Ne",
"Na", "Mg", "Al", "Si", "P" , "S" , "Cl", "Ar", "K" , "Ca",
"Sc", "Ti", "V" , "Cr", "Mn", "Fe", "Co", "Ni", "Cu", "Zn",
"Ga", "Ge", "As", "Se", "Br", "Kr", "Rb", "Sr", "Y" , "Zr",
"Nb", "Mo", "Tc", "Ru", "Rh", "Pd", "Ag", "Cd", "In", "Sn",
"Sb", "Te", "I" , "Xe", "Cs", "Ba", "La", "Ce", "Pr", "Nd",
"Pm", "Sm", "Eu", "Gd", "Tb", "Dy", "Ho", "Er", "Tm", "Yb",
"Lu", "Hf", "Ta", "W" , "Re", "Os", "Ir", "Pt", "Au", "Hg",
"Tl", "Pb", "Bi", "Po", "At", "Rn", "Fr", "Ra", "Ac", "Th",
"Pa", "U" , "Np", "Pu", "Am", "Cm", "Bk", "Cf", "Es", "Fm",
"Md", "No", "Lr", "Rf", "Db", "Sg", "Bh", "Hs", "Mt", "Ds",
"Rg", "Cn", "Nh", "Fl", "Mc", "Lv", "Ts", "Og",
"H", "He", "Li", "Be", "B", "C", "N", "O", "F", "Ne", "Na", "Mg", "Al", "Si", "P",
"S", "Cl", "Ar", "K", "Ca", "Sc", "Ti", "V", "Cr", "Mn", "Fe", "Co", "Ni", "Cu", "Zn",
"Ga", "Ge", "As", "Se", "Br", "Kr", "Rb", "Sr", "Y", "Zr", "Nb", "Mo", "Tc", "Ru", "Rh",
"Pd", "Ag", "Cd", "In", "Sn", "Sb", "Te", "I", "Xe", "Cs", "Ba", "La", "Ce", "Pr", "Nd",
"Pm", "Sm", "Eu", "Gd", "Tb", "Dy", "Ho", "Er", "Tm", "Yb", "Lu", "Hf", "Ta", "W", "Re",
"Os", "Ir", "Pt", "Au", "Hg", "Tl", "Pb", "Bi", "Po", "At", "Rn", "Fr", "Ra", "Ac", "Th",
"Pa", "U", "Np", "Pu", "Am", "Cm", "Bk", "Cf", "Es", "Fm", "Md", "No", "Lr", "Rf", "Db",
"Sg", "Bh", "Hs", "Mt", "Ds", "Rg", "Cn", "Nh", "Fl", "Mc", "Lv", "Ts", "Og",
};
int ntypes = atom->ntypes;
@ -121,8 +117,7 @@ FixMDIQM::FixMDIQM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
int anum;
for (anum = 0; anum < MAXELEMENT; anum++)
if (strcmp(arg[iarg + i], symbols[anum]) == 0) break;
if (anum == MAXELEMENT)
error->all(FLERR,"Invalid chemical element in fix mdi/qm command");
if (anum == MAXELEMENT) error->all(FLERR, "Invalid chemical element in fix mdi/qm command");
elements[i] = anum + 1;
}
iarg += ntypes + 1;
@ -210,8 +205,7 @@ FixMDIQM::FixMDIQM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
sumflag = 0;
int nlocal = atom->nlocal;
for (int i = 0; i < nlocal; i++)
array_atom[i][0] = array_atom[i][1] = array_atom[i][2] = 0.0;
for (int i = 0; i < nlocal; i++) array_atom[i][0] = array_atom[i][1] = array_atom[i][2] = 0.0;
}
/* ---------------------------------------------------------------------- */
@ -325,6 +319,14 @@ void FixMDIQM::init()
ierr = MDI_Check_command_exists("@DEFAULT", "<STRESS", mdicomm, &stress_exists);
if (ierr) error->all(FLERR, "MDI: <STRESS command check");
MPI_Bcast(&stress_exists, 1, MPI_INT, 0, world);
ierr = MDI_Check_command_exists("@DEFAULT", "<PE", mdicomm, &pe_exists);
if (ierr) error->all(FLERR, "MDI: <PE command check");
MPI_Bcast(&pe_exists, 1, MPI_INT, 0, world);
ierr = MDI_Check_command_exists("@DEFAULT", "<KE_ELEC", mdicomm, &keelec_exists);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC command check");
MPI_Bcast(&keelec_exists, 1, MPI_INT, 0, world);
}
// extract pointers to MC variables from id_mcfix
@ -337,7 +339,8 @@ void FixMDIQM::init()
int dim;
mc_active_ptr = (int *) f_mc->extract("mc_active", dim);
if (!mc_active_ptr || dim != 0)
error->all(FLERR,"Fix mdi/qm could not query mc_active from Monte Carlo fix ID {}", id_mcfix);
error->all(FLERR, "Fix mdi/qm could not query mc_active from Monte Carlo fix ID {}",
id_mcfix);
exclusion_group_ptr = (int *) f_mc->extract("exclusion_group", dim);
}
@ -363,7 +366,8 @@ void FixMDIQM::init()
// check if box has changed
if (new_system) set_box();
if (new_system)
set_box();
else {
double old_cell[9], old_cell_displ[3];
memcpy(old_cell, qm_cell, 9 * sizeof(double));
@ -378,7 +382,8 @@ void FixMDIQM::init()
// check if atom elements or types have changed
if (elements && elements_exists) {
if (new_system) set_eqm();
if (new_system)
set_eqm();
else {
int *eqm_old;
memory->create(eqm_old, nqm, "mdi/qm:eqm_old");
@ -390,7 +395,8 @@ void FixMDIQM::init()
}
} else if (types_exists) {
if (new_system) set_tqm();
if (new_system)
set_tqm();
else {
int *tqm_old;
memory->create(tqm_old, nqm, "mdi/qm:tqm_old");
@ -409,8 +415,10 @@ void FixMDIQM::init()
if (new_system) {
send_natoms();
send_box();
if (elements && elements_exists) send_elements();
else if (types_exists) send_types();
if (elements && elements_exists)
send_elements();
else if (types_exists)
send_types();
nqm_last = nqm;
}
}
@ -460,8 +468,10 @@ void FixMDIQM::post_force(int vflag)
// incremental when a system is slowly evolving (AIMD)
int new_system = 0;
if (nqm != nqm_last) new_system = 1;
else if (mcflag && *mc_active_ptr) new_system = 1;
if (nqm != nqm_last)
new_system = 1;
else if (mcflag && *mc_active_ptr)
new_system = 1;
// send new system info to MDI engine: atom count and elements/types
// reset QM data structures if atom count has changed
@ -502,16 +512,11 @@ void FixMDIQM::post_force(int vflag)
ierr = MDI_Send(&xqm[0][0], 3 * nqm, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: >COORDS data");
// request potential energy from MDI engine
// request QM energy from MDI engine
// this triggers engine to perform QM calculation
// qm_energy = fix output for global QM energy
// sets qm_energy = fix output for global QM energy
ierr = MDI_Send_command("<PE", mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE command");
ierr = MDI_Recv(&qm_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE data");
MPI_Bcast(&qm_energy, 1, MPI_DOUBLE, 0, world);
qm_energy *= mdi2lmp_energy;
request_qm_energy();
// request forces from MDI engine
@ -698,8 +703,7 @@ int FixMDIQM::set_nqm()
{
// require 3*nqm be a small INT, so can MPI_Allreduce xqm
if (3*atom->natoms > MAXSMALLINT)
error->all(FLERR,"Fix mdi/qm has too many atoms");
if (3 * atom->natoms > MAXSMALLINT) error->all(FLERR, "Fix mdi/qm has too many atoms");
int ncount = atom->natoms;
nexclude = 0;
@ -748,8 +752,10 @@ void FixMDIQM::create_qm_list()
int nqm_mine = 0;
for (int i = 0; i < nlocal; i++) {
if (!nexclude) nqm_mine++;
else if (!(mask[i] & excludebit)) nqm_mine++;
if (!nexclude)
nqm_mine++;
else if (!(mask[i] & excludebit))
nqm_mine++;
}
tagint *qmIDs_mine;
@ -757,8 +763,10 @@ void FixMDIQM::create_qm_list()
nqm_mine = 0;
for (int i = 0; i < nlocal; i++) {
if (!nexclude) qmIDs_mine[nqm_mine++] = tag[i];
else if (!(mask[i] & excludebit)) qmIDs_mine[nqm_mine++] = tag[i];
if (!nexclude)
qmIDs_mine[nqm_mine++] = tag[i];
else if (!(mask[i] & excludebit))
qmIDs_mine[nqm_mine++] = tag[i];
}
int *recvcounts, *displs;
@ -771,8 +779,8 @@ void FixMDIQM::create_qm_list()
for (int iproc = 1; iproc < nprocs; iproc++)
displs[iproc] = displs[iproc - 1] + recvcounts[iproc - 1];
MPI_Allgatherv(qmIDs_mine,nqm_mine,MPI_LMP_TAGINT,qmIDs,recvcounts,displs,
MPI_LMP_TAGINT,world);
MPI_Allgatherv(qmIDs_mine, nqm_mine, MPI_LMP_TAGINT, qmIDs, recvcounts, displs, MPI_LMP_TAGINT,
world);
memory->destroy(qmIDs_mine);
memory->destroy(recvcounts);
@ -818,8 +826,10 @@ void FixMDIQM::set_qm2owned()
for (int i = 0; i < nqm; i++) {
index = atom->map(qmIDs[i]);
if (index >= nlocal) qm2owned[i] = -1;
else qm2owned[i] = index;
if (index >= nlocal)
qm2owned[i] = -1;
else
qm2owned[i] = index;
}
}
@ -974,6 +984,9 @@ void FixMDIQM::send_box()
{
int ierr;
// only send cell dimensions if fully periodic simulation
if (domain->nonperiodic == 0) {
if (celldispl_exists) {
ierr = MDI_Send_command(">CELL_DISPL", mdicomm);
if (ierr) error->all(FLERR, "MDI: >CELL_DISPL command");
@ -985,6 +998,48 @@ void FixMDIQM::send_box()
if (ierr) error->all(FLERR, "MDI: >CELL command");
ierr = MDI_Send(qm_cell, 9, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: >CELL data");
} else if (domain->xperiodic == 1 || domain->yperiodic == 1 ||
domain->zperiodic == 1) {
error->all(FLERR,"MDI requires fully periodic or fully non-periodic system");
}
}
/* ----------------------------------------------------------------------
request QM energy from MDI engine
set qm_energy = fix output for global QM energy
------------------------------------------------------------------------- */
void FixMDIQM::request_qm_energy()
{
int ierr;
// QM energy = <PE + <KE_ELEC or <ENERGY, depending on engine options
if (pe_exists && keelec_exists) {
int pe_energy, keelec_energy;
ierr = MDI_Send_command("<PE", mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE command");
ierr = MDI_Recv(&pe_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE data");
ierr = MDI_Send_command("<KE_ELEC", mdicomm);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC command");
ierr = MDI_Recv(&keelec_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC data");
qm_energy = pe_energy + keelec_energy;
} else {
ierr = MDI_Send_command("<ENERGY", mdicomm);
if (ierr) error->all(FLERR, "MDI: <ENERGY command");
ierr = MDI_Recv(&qm_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <ENERGY data");
}
MPI_Bcast(&qm_energy, 1, MPI_DOUBLE, 0, world);
qm_energy *= mdi2lmp_energy;
}
/* ----------------------------------------------------------------------

View File

@ -62,7 +62,7 @@ class FixMDIQM : public Fix {
MDI_Comm mdicomm;
int natoms_exists, celldispl_exists, elements_exists, types_exists;
int stress_exists;
int stress_exists, pe_exists, keelec_exists;
int nmax;
@ -105,6 +105,8 @@ class FixMDIQM : public Fix {
void send_elements();
void send_box();
void request_qm_energy();
void unit_conversions();
};

View File

@ -44,8 +44,7 @@ FixMDIQMMM::FixMDIQMMM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
// check requirements for LAMMPS to work with MDI for a QMMM engine
// atom IDs do not need to be consecutive
if (!atom->tag_enable)
error->all(FLERR,"Fix mdi/qmmm requires atom IDs be defined");
if (!atom->tag_enable) error->all(FLERR, "Fix mdi/qmmm requires atom IDs be defined");
if (atom->map_style == Atom::MAP_NONE)
error->all(FLERR, "Fix mdi/qmmm requires an atom map be defined");
@ -59,9 +58,12 @@ FixMDIQMMM::FixMDIQMMM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
// mode arg
if (strcmp(arg[3],"direct") == 0) mode = DIRECT;
else if (strcmp(arg[3],"potential") == 0) mode = POTENTIAL;
else error->all(FLERR,"Illegal fix mdi/qmmm command");
if (strcmp(arg[3], "direct") == 0)
mode = DIRECT;
else if (strcmp(arg[3], "potential") == 0)
mode = POTENTIAL;
else
error->all(FLERR, "Illegal fix mdi/qmmm command");
// optional args
@ -92,18 +94,14 @@ FixMDIQMMM::FixMDIQMMM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
} else if (strcmp(arg[iarg], "elements") == 0) {
const char *symbols[] = {
"H" , "He", "Li", "Be", "B" , "C" , "N" , "O" , "F" , "Ne",
"Na", "Mg", "Al", "Si", "P" , "S" , "Cl", "Ar", "K" , "Ca",
"Sc", "Ti", "V" , "Cr", "Mn", "Fe", "Co", "Ni", "Cu", "Zn",
"Ga", "Ge", "As", "Se", "Br", "Kr", "Rb", "Sr", "Y" , "Zr",
"Nb", "Mo", "Tc", "Ru", "Rh", "Pd", "Ag", "Cd", "In", "Sn",
"Sb", "Te", "I" , "Xe", "Cs", "Ba", "La", "Ce", "Pr", "Nd",
"Pm", "Sm", "Eu", "Gd", "Tb", "Dy", "Ho", "Er", "Tm", "Yb",
"Lu", "Hf", "Ta", "W" , "Re", "Os", "Ir", "Pt", "Au", "Hg",
"Tl", "Pb", "Bi", "Po", "At", "Rn", "Fr", "Ra", "Ac", "Th",
"Pa", "U" , "Np", "Pu", "Am", "Cm", "Bk", "Cf", "Es", "Fm",
"Md", "No", "Lr", "Rf", "Db", "Sg", "Bh", "Hs", "Mt", "Ds",
"Rg", "Cn", "Nh", "Fl", "Mc", "Lv", "Ts", "Og",
"H", "He", "Li", "Be", "B", "C", "N", "O", "F", "Ne", "Na", "Mg", "Al", "Si", "P",
"S", "Cl", "Ar", "K", "Ca", "Sc", "Ti", "V", "Cr", "Mn", "Fe", "Co", "Ni", "Cu", "Zn",
"Ga", "Ge", "As", "Se", "Br", "Kr", "Rb", "Sr", "Y", "Zr", "Nb", "Mo", "Tc", "Ru", "Rh",
"Pd", "Ag", "Cd", "In", "Sn", "Sb", "Te", "I", "Xe", "Cs", "Ba", "La", "Ce", "Pr", "Nd",
"Pm", "Sm", "Eu", "Gd", "Tb", "Dy", "Ho", "Er", "Tm", "Yb", "Lu", "Hf", "Ta", "W", "Re",
"Os", "Ir", "Pt", "Au", "Hg", "Tl", "Pb", "Bi", "Po", "At", "Rn", "Fr", "Ra", "Ac", "Th",
"Pa", "U", "Np", "Pu", "Am", "Cm", "Bk", "Cf", "Es", "Fm", "Md", "No", "Lr", "Rf", "Db",
"Sg", "Bh", "Hs", "Mt", "Ds", "Rg", "Cn", "Nh", "Fl", "Mc", "Lv", "Ts", "Og",
};
int ntypes = atom->ntypes;
@ -222,8 +220,7 @@ FixMDIQMMM::FixMDIQMMM(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg)
sumflag = 0;
int nlocal = atom->nlocal;
for (int i = 0; i < nlocal; i++)
array_atom[i][0] = array_atom[i][1] = array_atom[i][2] = 0.0;
for (int i = 0; i < nlocal; i++) array_atom[i][0] = array_atom[i][1] = array_atom[i][2] = 0.0;
}
/* ---------------------------------------------------------------------- */
@ -360,6 +357,14 @@ void FixMDIQMMM::init()
if (ierr) error->all(FLERR, "MDI: <STRESS command check");
MPI_Bcast(&stress_exists, 1, MPI_INT, 0, world);
ierr = MDI_Check_command_exists("@DEFAULT", "<PE", mdicomm, &pe_exists);
if (ierr) error->all(FLERR, "MDI: <PE command check");
MPI_Bcast(&pe_exists, 1, MPI_INT, 0, world);
ierr = MDI_Check_command_exists("@DEFAULT", "<KE_ELEC", mdicomm, &keelec_exists);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC command check");
MPI_Bcast(&keelec_exists, 1, MPI_INT, 0, world);
// check if MDI engine supports commands that match QMMM mode
if (mode == POTENTIAL) {
@ -425,7 +430,8 @@ void FixMDIQMMM::init()
pair_coul = force->pair_match("coul/cut", 1, 0);
if (!pair_coul) pair_coul = force->pair_match("coul/long", 1, 0);
if (!pair_coul) pair_coul = force->pair_match("coul/msm", 1, 0);
if (!pair_coul) error->all(FLERR,"Fix mdi/qmmm potential requires Coulomb-only pair sub-style");
if (!pair_coul)
error->all(FLERR, "Fix mdi/qmmm potential requires Coulomb-only pair sub-style");
}
// determine whether a new vs incremental QMMM calc is needed
@ -460,7 +466,8 @@ void FixMDIQMMM::init()
// check if box has changed
if (new_system) set_box();
if (new_system)
set_box();
else {
double old_cell[9], old_cell_displ[3];
memcpy(old_cell, qm_cell, 9 * sizeof(double));
@ -530,14 +537,18 @@ void FixMDIQMMM::init()
if (new_system) {
send_natoms_qm();
send_box();
if (elements && elements_exists) send_elements_qm();
else if (types_exists) send_types_qm();
if (elements && elements_exists)
send_elements_qm();
else if (types_exists)
send_types_qm();
nqm_last = nqm;
if (mode == DIRECT) {
send_natoms_mm();
if (elements && elements_exists) send_elements_mm();
else if (types_exists) send_types_mm();
if (elements && elements_exists)
send_elements_mm();
else if (types_exists)
send_types_mm();
set_qmm();
send_charges_mm();
}
@ -614,14 +625,12 @@ void FixMDIQMMM::pre_force(int vflag)
int ntotal = nlocal;
if (force->newton_pair) ntotal += atom->nghost;
for (int i = 0; i < ntotal; i++)
ecoul[i] = eatom_pair[i];
for (int i = 0; i < ntotal; i++) ecoul[i] = eatom_pair[i];
if (force->newton_pair) comm->reverse_comm(this);
if (force->kspace) {
for (int i = 0; i < nlocal; i++)
ecoul[i] += eatom_kspace[i];
for (int i = 0; i < nlocal; i++) ecoul[i] += eatom_kspace[i];
}
// setup QM inputs: xqm and qpotential
@ -644,8 +653,10 @@ void FixMDIQMMM::pre_force(int vflag)
for (int i = 0; i < nqm; i++) {
ilocal = qm2owned[i];
if (ilocal >= 0) {
if (q[ilocal] == 0.0) qpotential_mine[i] = 0.0;
else qpotential_mine[i] = 2.0 * ecoul[ilocal] / q[ilocal];
if (q[ilocal] == 0.0)
qpotential_mine[i] = 0.0;
else
qpotential_mine[i] = 2.0 * ecoul[ilocal] / q[ilocal];
}
}
@ -715,16 +726,11 @@ void FixMDIQMMM::pre_force(int vflag)
ierr = MDI_Send(qpotential, nqm, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: >POTENTIAL_AT_NUCLEI data");
// request QM potential energy from MDI engine
// request QM energy from MDI engine
// this triggers engine to perform QM calculation
// qm_energy = fix output for global QM energy
// sets qm_energy = fix output for global QM energy
ierr = MDI_Send_command("<PE", mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE command");
ierr = MDI_Recv(&qm_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE data");
MPI_Bcast(&qm_energy, 1, MPI_DOUBLE, 0, world);
qm_energy *= mdi2lmp_energy;
request_qm_energy();
// request forces on QM atoms from MDI engine
@ -826,8 +832,10 @@ void FixMDIQMMM::pre_force(int vflag)
void FixMDIQMMM::post_force(int vflag)
{
if (mode == DIRECT) post_force_direct(vflag);
else if (mode == POTENTIAL) post_force_potential(vflag);
if (mode == DIRECT)
post_force_direct(vflag);
else if (mode == POTENTIAL)
post_force_potential(vflag);
}
/* ----------------------------------------------------------------------
@ -878,16 +886,11 @@ void FixMDIQMMM::post_force_direct(int vflag)
ierr = MDI_Send(&xmm[0][0], 3 * nmm, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: >CLATTICE data");
// request QM potential energy from MDI engine
// request QM energy from MDI engine
// this triggers engine to perform QM calculation
// qm_energy = fix output for global QM energy
// sets qm_energy = fix output for global QM energy
ierr = MDI_Send_command("<PE", mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE command");
ierr = MDI_Recv(&qm_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE data");
MPI_Bcast(&qm_energy, 1, MPI_DOUBLE, 0, world);
qm_energy *= mdi2lmp_energy;
request_qm_energy();
// request forces on QM atoms from MDI engine
@ -1123,8 +1126,7 @@ void FixMDIQMMM::min_post_force(int vflag)
/* ---------------------------------------------------------------------- */
int FixMDIQMMM::pack_forward_comm(int n, int *list, double *buf,
int /*pbc_flag*/, int * /*pbc*/)
int FixMDIQMMM::pack_forward_comm(int n, int *list, double *buf, int /*pbc_flag*/, int * /*pbc*/)
{
int i, j, m;
@ -1304,8 +1306,7 @@ int FixMDIQMMM::set_nqm()
// require 3*nqm be a small INT, so can MPI_Allreduce QM values
if (3*ngroup > MAXSMALLINT)
error->all(FLERR,"Fix mdi/qmmm has too many quantum atoms");
if (3 * ngroup > MAXSMALLINT) error->all(FLERR, "Fix mdi/qmmm has too many quantum atoms");
// error if nqm = 0
// error if nqm = natoms, should use fix mdi/qm instead
@ -1369,8 +1370,8 @@ void FixMDIQMMM::create_qm_list()
for (int iproc = 1; iproc < nprocs; iproc++)
displs[iproc] = displs[iproc - 1] + recvcounts[iproc - 1];
MPI_Allgatherv(qmIDs_mine,nqm_mine,MPI_LMP_TAGINT,qmIDs,recvcounts,displs,
MPI_LMP_TAGINT,world);
MPI_Allgatherv(qmIDs_mine, nqm_mine, MPI_LMP_TAGINT, qmIDs, recvcounts, displs, MPI_LMP_TAGINT,
world);
memory->destroy(qmIDs_mine);
memory->destroy(recvcounts);
@ -1436,8 +1437,8 @@ void FixMDIQMMM::create_mm_list()
for (int iproc = 1; iproc < nprocs; iproc++)
displs[iproc] = displs[iproc - 1] + recvcounts[iproc - 1];
MPI_Allgatherv(mmIDs_mine,nmm_mine,MPI_LMP_TAGINT,mmIDs,recvcounts,displs,
MPI_LMP_TAGINT,world);
MPI_Allgatherv(mmIDs_mine, nmm_mine, MPI_LMP_TAGINT, mmIDs, recvcounts, displs, MPI_LMP_TAGINT,
world);
memory->destroy(mmIDs_mine);
memory->destroy(recvcounts);
@ -1483,8 +1484,10 @@ void FixMDIQMMM::set_qm2owned()
for (int i = 0; i < nqm; i++) {
index = atom->map(qmIDs[i]);
if (index >= nlocal) qm2owned[i] = -1;
else qm2owned[i] = index;
if (index >= nlocal)
qm2owned[i] = -1;
else
qm2owned[i] = index;
}
}
@ -1501,8 +1504,10 @@ void FixMDIQMMM::set_mm2owned()
for (int i = 0; i < nmm; i++) {
index = atom->map(mmIDs[i]);
if (index >= nlocal) mm2owned[i] = -1;
else mm2owned[i] = index;
if (index >= nlocal)
mm2owned[i] = -1;
else
mm2owned[i] = index;
}
}
@ -1529,7 +1534,6 @@ void FixMDIQMMM::set_box()
for (int icell = 0; icell < 9; icell++) qm_cell[icell] *= lmp2mdi_length;
}
/* ----------------------------------------------------------------------
fill xqm with QM atom coords
if convert, perform LAMMPS to MDI unit conversion
@ -1717,7 +1721,6 @@ void FixMDIQMMM::send_natoms_qm()
ierr = MDI_Send(&nqm, 1, MDI_INT, mdicomm);
if (ierr) error->all(FLERR, "MDI: >NATOMS data");
} else {
ierr = MDI_Send_command("<NATOMS", mdicomm);
if (ierr) error->all(FLERR, "MDI: <NATOMS command");
@ -1764,6 +1767,9 @@ void FixMDIQMMM::send_box()
{
int ierr;
// only send cell dimensions if fully periodic simulation
if (domain->nonperiodic == 0) {
if (celldispl_exists) {
ierr = MDI_Send_command(">CELL_DISPL", mdicomm);
if (ierr) error->all(FLERR, "MDI: >CELL_DISPL command");
@ -1775,6 +1781,11 @@ void FixMDIQMMM::send_box()
if (ierr) error->all(FLERR, "MDI: >CELL command");
ierr = MDI_Send(qm_cell, 9, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: >CELL data");
} else if (domain->xperiodic == 1 || domain->yperiodic == 1 ||
domain->zperiodic == 1) {
error->all(FLERR,"MDI requires fully periodic or fully non-periodic system");
}
}
/* ----------------------------------------------------------------------
@ -1825,6 +1836,43 @@ void FixMDIQMMM::send_charges_mm()
if (ierr) error->all(FLERR, "MDI: >LATTICE data");
}
/* ----------------------------------------------------------------------
request QM energy from MDI engine
set qm_energy = fix output for global QM energy
------------------------------------------------------------------------- */
void FixMDIQMMM::request_qm_energy()
{
int ierr;
// QM energy = <PE + <KE_ELEC or <ENERGY, depending on engine options
if (pe_exists && keelec_exists) {
int pe_energy, keelec_energy;
ierr = MDI_Send_command("<PE", mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE command");
ierr = MDI_Recv(&pe_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <PE data");
ierr = MDI_Send_command("<KE_ELEC", mdicomm);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC command");
ierr = MDI_Recv(&keelec_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC data");
qm_energy = pe_energy + keelec_energy;
} else {
ierr = MDI_Send_command("<ENERGY", mdicomm);
if (ierr) error->all(FLERR, "MDI: <ENERGY command");
ierr = MDI_Recv(&qm_energy, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <ENERGY data");
}
MPI_Bcast(&qm_energy, 1, MPI_DOUBLE, 0, world);
qm_energy *= mdi2lmp_energy;
}
/* ----------------------------------------------------------------------
MDI to/from LAMMPS conversion factors
------------------------------------------------------------------------- */

View File

@ -68,7 +68,7 @@ class FixMDIQMMM : public Fix {
MDI_Comm mdicomm;
int natoms_exists, celldispl_exists, elements_exists, types_exists;
int stress_exists;
int stress_exists, pe_exists, keelec_exists;
int nmax;
@ -156,8 +156,9 @@ class FixMDIQMMM : public Fix {
void send_elements_mm();
void send_charges_mm();
void unit_conversions();
void request_qm_energy();
void unit_conversions();
};
} // namespace LAMMPS_NS

View File

@ -430,6 +430,9 @@ int MDIEngine::execute_command(const char *command, MDI_Comm mdicomm)
if (!actionflag && strcmp(node_engine, "@DEFAULT") == 0) evaluate();
send_pe();
} else if (strcmp(command, "<KE_ELEC") == 0) {
send_ke_elec();
} else if (strcmp(command, "<STRESS") == 0) {
if (!actionflag && strcmp(node_engine, "@DEFAULT") == 0) evaluate();
send_stress();
@ -540,6 +543,7 @@ void MDIEngine::mdi_commands()
MDI_Register_command("@DEFAULT", "<MASSES");
MDI_Register_command("@DEFAULT", "<NATOMS");
MDI_Register_command("@DEFAULT", "<PE");
MDI_Register_command("@DEFAULT", "<KE_ELEC");
MDI_Register_command("@DEFAULT", "<STRESS");
MDI_Register_command("@DEFAULT", "<TYPES");
MDI_Register_command("@DEFAULT", "<VELOCITIES");
@ -618,6 +622,7 @@ void MDIEngine::mdi_commands()
MDI_Register_command("@FORCES", "<FORCES");
MDI_Register_command("@FORCES", "<KE");
MDI_Register_command("@FORCES", "<PE");
MDI_Register_command("@FORCES", "<KE_ELEC");
MDI_Register_command("@FORCES", "<STRESS");
MDI_Register_command("@FORCES", "<VELOCITIES");
MDI_Register_command("@FORCES", ">FORCES");
@ -639,6 +644,7 @@ void MDIEngine::mdi_commands()
MDI_Register_command("@ENDSTEP", "<FORCES");
MDI_Register_command("@ENDSTEP", "<KE");
MDI_Register_command("@ENDSTEP", "<PE");
MDI_Register_command("@ENDSTEP", "<KE_ELEC");
MDI_Register_command("@ENDSTEP", "<STRESS");
MDI_Register_command("@ENDSTEP", "@");
MDI_Register_command("@ENDSTEP", "@DEFAULT");
@ -1525,6 +1531,21 @@ void MDIEngine::send_pe()
if (ierr) error->all(FLERR, "MDI: <PE data");
}
/* ----------------------------------------------------------------------
<KE_ELEC command
send kinetic energy of the electrons
zero for LAMMPS, because it does not model electrons explicitly
for compatibiity with QM engines which support this command
---------------------------------------------------------------------- */
void MDIEngine::send_ke_elec()
{
double ke_elec = 0.0;
int ierr = MDI_Send(&ke_elec, 1, MDI_DOUBLE, mdicomm);
if (ierr) error->all(FLERR, "MDI: <KE_ELEC data");
}
/* ----------------------------------------------------------------------
<STRESS command
send 9-component stress tensor (no kinetic energy term)

View File

@ -125,6 +125,7 @@ class MDIEngine : protected Pointers {
void send_labels();
void send_natoms();
void send_pe();
void send_ke_elec();
void send_stress();
void send_double1(int);