#include <mpi.h>
#include <cmath>
#include <cstring>
#include "dynamical_matrix.h"
#include "atom.h"
#include "domain.h"
#include "comm.h"
#include "error.h"
#include "group.h"
#include "force.h"
#include "memory.h"
#include "bond.h"
#include "angle.h"
#include "dihedral.h"
#include "improper.h"
#include "kspace.h"
#include "update.h"
#include "neighbor.h"
#include "pair.h"
#include "timer.h"
#include "finish.h"
#include <algorithm>
using namespace LAMMPS_NS;
enum{REGULAR,ESKM};
DynamicalMatrix::DynamicalMatrix(LAMMPS *lmp) : Pointers(lmp), fp(NULL)
{
external_force_clear = 1;
}
DynamicalMatrix::~DynamicalMatrix()
{
if (fp && me == 0) fclose(fp);
memory->destroy(groupmap);
fp = NULL;
}
void DynamicalMatrix::setup()
{
if (triclinic) domain->x2lamda(atom->nlocal);
domain->pbc();
domain->reset_box();
comm->setup();
if (neighbor->style) neighbor->setup_bins();
comm->exchange();
comm->borders();
if (triclinic) domain->lamda2x(atom->nlocal+atom->nghost);
domain->image_check();
domain->box_too_small_check();
neighbor->build(1);
neighbor->ncalls = 0;
neighbor->every = 2; neighbor->delay = 1;
neighbor->ago = 0;
neighbor->ndanger = 0;
external_force_clear = 0;
eflag=0;
vflag=0;
update_force();
if (gcount == atom->natoms)
for (bigint i=0; i<atom->natoms; i++)
groupmap[i] = i;
else
create_groupmap();
}
void DynamicalMatrix::command(int narg, char **arg)
{
MPI_Comm_rank(world,&me);
if (domain->box_exist == 0)
error->all(FLERR,"Dynamical_matrix command before simulation box is defined");
if (narg < 2) error->all(FLERR,"Illegal dynamical_matrix command");
lmp->init();
triclinic = domain->triclinic;
if (force->pair && force->pair->compute_flag) pair_compute_flag = 1;
else pair_compute_flag = 0;
if (force->kspace && force->kspace->compute_flag) kspace_compute_flag = 1;
else kspace_compute_flag = 0;
igroup = group->find(arg[0]);
if (igroup == -1) error->all(FLERR,"Could not find dynamical matrix group ID");
groupbit = group->bitmask[igroup];
gcount = group->count(igroup);
dynlen = (gcount)*3;
memory->create(groupmap,atom->natoms,"total_group_map:totalgm");
update->setupflag = 1;
int style = -1;
if (strcmp(arg[1],"regular") == 0) style = REGULAR;
else if (strcmp(arg[1],"eskm") == 0) style = ESKM;
else error->all(FLERR,"Illegal Dynamical Matrix command");
del = force->numeric(FLERR, arg[2]);
binaryflag = 0;
scaleflag = 0;
compressed = 0;
file_flag = 0;
file_opened = 0;
conversion = 1;
if (style == REGULAR) options(narg-3,&arg[3]); else if (style == ESKM) options(narg-3,&arg[3]); else if (comm->me == 0 && screen) fprintf(screen,"Illegal Dynamical Matrix command\n");
if (atom->map_style == 0)
error->all(FLERR,"Dynamical_matrix command requires an atom map, see atom_modify");
if (style == REGULAR) {
setup();
timer->init();
timer->barrier_start();
calculateMatrix();
timer->barrier_stop();
}
if (style == ESKM) {
setup();
convert_units(update->unit_style);
conversion = conv_energy/conv_distance/conv_mass;
timer->init();
timer->barrier_start();
calculateMatrix();
timer->barrier_stop();
}
Finish finish(lmp);
finish.end(1);
}
void DynamicalMatrix::options(int narg, char **arg)
{
if (narg < 0) error->all(FLERR,"Illegal dynamical_matrix command");
int iarg = 0;
const char* filename = "dynmat.dyn";
while (iarg < narg) {
if (strcmp(arg[iarg],"binary") == 0) {
if (iarg + 2 > narg) error->all(FLERR, "Illegal dynamical_matrix command");
if (strcmp(arg[iarg+1],"gzip") == 0) {
compressed = 1;
}
else if (strcmp(arg[iarg+1],"yes") == 0) {
binaryflag = 1;
}
iarg += 2;
}
else if (strcmp(arg[iarg],"file") == 0) {
if (iarg+2 > narg) error->all(FLERR, "Illegal dynamical_matrix command");
filename = arg[iarg + 1];
file_flag = 1;
iarg += 2;
} else error->all(FLERR,"Illegal dynamical_matrix command");
}
if (file_flag == 1) {
openfile(filename);
}
}
void DynamicalMatrix::openfile(const char* filename)
{
if (file_opened) return;
if (compressed) {
#ifdef LAMMPS_GZIP
char gzip[128];
sprintf(gzip,"gzip -6 > %s",filename);
#ifdef _WIN32
fp = _popen(gzip,"wb");
#else
fp = popen(gzip,"w");
#endif
#else
error->one(FLERR,"Cannot open gzipped file");
#endif
} else if (binaryflag) {
fp = fopen(filename,"wb");
} else {
fp = fopen(filename,"w");
}
if (fp == NULL) error->one(FLERR,"Cannot open dump file");
file_opened = 1;
}
void DynamicalMatrix::calculateMatrix()
{
int local_idx; int local_jdx; int nlocal = atom->nlocal;
bigint natoms = atom->natoms;
int *type = atom->type;
bigint *gm = groupmap;
double imass; double *m = atom->mass;
double **f = atom->f;
double **dynmat = new double*[3];
for (int i=0; i<3; i++)
dynmat[i] = new double[dynlen];
double **fdynmat = new double*[3];
for (int i=0; i<3; i++)
fdynmat[i] = new double[dynlen];
dynmat_clear(dynmat);
if (comm->me == 0 && screen) {
fprintf(screen,"Calculating Dynamical Matrix ...\n");
fprintf(screen," Total # of atoms = " BIGINT_FORMAT "\n", natoms);
fprintf(screen," Atoms in group = " BIGINT_FORMAT "\n", gcount);
fprintf(screen," Total dynamical matrix elements = " BIGINT_FORMAT "\n", (dynlen*dynlen) );
}
update->nsteps = 0;
int prog = 0;
for (bigint i=1; i<=natoms; i++){
local_idx = atom->map(i);
if (gm[i-1] < 0)
continue;
for (bigint alpha=0; alpha<3; alpha++){
displace_atom(local_idx, alpha, 1);
update_force();
for (bigint j=1; j<=natoms; j++){
local_jdx = atom->map(j);
if (local_idx >= 0 && local_jdx >= 0 && local_jdx < nlocal
&& gm[j-1] >= 0){
for (int beta=0; beta<3; beta++){
dynmat[alpha][gm[j-1]*3+beta] -= f[local_jdx][beta];
}
}
}
displace_atom(local_idx,alpha,-2);
update_force();
for (bigint j=1; j<=natoms; j++){
local_jdx = atom->map(j);
if (local_idx >= 0 && local_jdx >= 0 && local_jdx < nlocal
&& gm[j-1] >= 0){
for (bigint beta=0; beta<3; beta++){
if (atom->rmass_flag == 1)
imass = sqrt(m[local_idx] * m[local_jdx]);
else
imass = sqrt(m[type[local_idx]] * m[type[local_jdx]]);
dynmat[alpha][gm[j-1]*3+beta] -= -f[local_jdx][beta];
dynmat[alpha][gm[j-1]*3+beta] /= (2 * del * imass);
dynmat[alpha][gm[j-1]*3+beta] *= conversion;
}
}
}
displace_atom(local_idx,alpha,1);
}
for (int k=0; k<3; k++)
MPI_Reduce(dynmat[k],fdynmat[k],dynlen,MPI_DOUBLE,MPI_SUM,0,world);
if (me == 0)
writeMatrix(fdynmat);
dynmat_clear(dynmat);
if (comm->me == 0 && screen) {
int p = 10 * gm[i-1] / gcount;
if (p > prog) {
prog = p;
fprintf(screen," %d%%",p*10);
fflush(screen);
}
}
}
if (comm->me == 0 && screen) fprintf(screen,"\n");
for (int i=0; i < 3; i++)
delete [] dynmat[i];
delete [] dynmat;
for (int i=0; i < 3; i++)
delete [] fdynmat[i];
delete [] fdynmat;
if (screen && me ==0 ) fprintf(screen,"Finished Calculating Dynamical Matrix\n");
}
void DynamicalMatrix::writeMatrix(double **dynmat)
{
if (me != 0 || !fp)
return;
clearerr(fp);
if (binaryflag) {
for (int i=0; i<3; i++)
fwrite(dynmat[i], sizeof(double), dynlen, fp);
if (ferror(fp))
error->one(FLERR, "Error writing to binary file");
} else {
for (int i = 0; i < 3; i++) {
for (bigint j = 0; j < dynlen; j++) {
if ((j+1)%3==0) fprintf(fp, "%4.8f\n", dynmat[i][j]);
else fprintf(fp, "%4.8f ",dynmat[i][j]);
}
}
if (ferror(fp))
error->one(FLERR,"Error writing to file");
}
}
void DynamicalMatrix::displace_atom(int local_idx, int direction, int magnitude)
{
if (local_idx < 0) return;
double **x = atom->x;
int *sametag = atom->sametag;
int j = local_idx;
x[local_idx][direction] += del*magnitude;
while (sametag[j] >= 0){
j = sametag[j];
x[j][direction] += del*magnitude;
}
}
void DynamicalMatrix::update_force()
{
force_clear();
if (pair_compute_flag) {
force->pair->compute(eflag,vflag);
timer->stamp(Timer::PAIR);
}
if (atom->molecular) {
if (force->bond) force->bond->compute(eflag,vflag);
if (force->angle) force->angle->compute(eflag,vflag);
if (force->dihedral) force->dihedral->compute(eflag,vflag);
if (force->improper) force->improper->compute(eflag,vflag);
timer->stamp(Timer::BOND);
}
if (kspace_compute_flag) {
force->kspace->compute(eflag,vflag);
timer->stamp(Timer::KSPACE);
}
if (force->newton) {
comm->reverse_comm();
timer->stamp(Timer::COMM);
}
++ update->nsteps;
}
void DynamicalMatrix::force_clear()
{
if (external_force_clear) return;
size_t nbytes = sizeof(double) * atom->nlocal;
if (force->newton) nbytes += sizeof(double) * atom->nghost;
if (nbytes) {
memset(&atom->f[0][0],0,3*nbytes);
}
}
void DynamicalMatrix::dynmat_clear(double **dynmat)
{
size_t nbytes = sizeof(double) * dynlen;
if (nbytes) {
for (int i=0; i<3; i++)
memset(&dynmat[i][0],0,nbytes);
}
}
void DynamicalMatrix::convert_units(const char *style)
{
if (strcmp(style,"lj") == 0) {
error->all(FLERR,"Conversion Not Set");
} else if (strcmp(style,"real") == 0) {
conv_energy = 418.4; conv_mass = 1; conv_distance = 1;
} else if (strcmp(style,"metal") == 0) {
conv_energy = 9648.5; conv_mass = 1; conv_distance = 1;
} else if (strcmp(style,"si") == 0) {
if (comm->me) error->warning(FLERR,"Conversion Warning: Multiplication by Large Float");
conv_energy = 6.022E22; conv_mass = 6.022E26; conv_distance = 1E-10;
} else if (strcmp(style,"cgs") == 0) {
if (comm->me) error->warning(FLERR,"Conversion Warning: Multiplication by Large Float");
conv_energy = 6.022E12; conv_mass = 6.022E23; conv_distance = 1E-7;
} else if (strcmp(style,"electron") == 0) {
conv_energy = 262550; conv_mass = 1; conv_distance = 0.529177249;
} else if (strcmp(style,"micro") == 0) {
if (comm->me) error->warning(FLERR,"Conversion Warning: Untested Conversion");
conv_energy = 6.022E10; conv_mass = 6.022E11; conv_distance = 1E-4;
} else if (strcmp(style,"nano") == 0) {
if (comm->me) error->warning(FLERR,"Conversion Warning: Untested Conversion");
conv_energy = 6.022E4; conv_mass = 6.022E5; conv_distance = 0.1;
} else error->all(FLERR,"Units Type Conversion Not Found");
}
void DynamicalMatrix::create_groupmap()
{
int local_idx; int gid = 0; int nlocal = atom->nlocal;
int *mask = atom->mask;
bigint natoms = atom->natoms;
int *recv = new int[comm->nprocs];
int *displs = new int[comm->nprocs];
bigint *temp_groupmap = new bigint[natoms];
for (bigint i=1; i<=natoms; i++){
local_idx = atom->map(i);
if ((local_idx >= 0) && (local_idx < nlocal) && mask[local_idx] & groupbit)
gid += 1; }
bigint *sub_groupmap = new bigint[gid];
gid = 0;
for (bigint i=1; i<=natoms; i++){
local_idx = atom->map(i);
if ((local_idx >= 0) && (local_idx < nlocal) && mask[local_idx] & groupbit){
sub_groupmap[gid] = i;
gid += 1;
}
}
for (int i=0; i<comm->nprocs; i++){
recv[i] = 0;
}
recv[comm->me] = gid;
MPI_Allreduce(recv,displs,comm->nprocs,MPI_INT,MPI_SUM,world);
for (int i=0; i<comm->nprocs; i++){
recv[i]=displs[i];
if (i>0) displs[i] = displs[i-1]+recv[i-1];
else displs[i] = 0;
}
MPI_Allgatherv(sub_groupmap,gid,MPI_LMP_BIGINT,temp_groupmap,recv,displs,MPI_LMP_BIGINT,world);
std::sort(temp_groupmap,temp_groupmap+gcount);
bigint j = 0;
for (bigint i=1; i<=natoms; i++){
if (j < gcount && i == temp_groupmap[j])
groupmap[i-1] = j++;
else
groupmap[i-1] = -1;
}
delete[] recv;
delete[] displs;
delete[] sub_groupmap;
delete[] temp_groupmap;
}