#include "pair_smd_triangulated_surface.h"
#include <mpi.h>
#include <cmath>
#include <cstdlib>
#include <cstring>
#include <Eigen/Eigen>
#include "atom.h"
#include "domain.h"
#include "force.h"
#include "comm.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "neigh_request.h"
#include "memory.h"
#include "error.h"
using namespace std;
using namespace LAMMPS_NS;
using namespace Eigen;
#define SQRT2 1.414213562e0
PairTriSurf::PairTriSurf(LAMMPS *lmp) :
Pair(lmp) {
onerad_dynamic = onerad_frozen = maxrad_dynamic = maxrad_frozen = NULL;
bulkmodulus = NULL;
kn = NULL;
scale = 1.0;
}
PairTriSurf::~PairTriSurf() {
if (allocated) {
memory->destroy(setflag);
memory->destroy(cutsq);
memory->destroy(bulkmodulus);
memory->destroy(kn);
delete[] onerad_dynamic;
delete[] onerad_frozen;
delete[] maxrad_dynamic;
delete[] maxrad_frozen;
}
}
void PairTriSurf::compute(int eflag, int vflag) {
int i, j, ii, jj, inum, jnum, itype, jtype;
double rsq, r, evdwl, fpair;
int *ilist, *jlist, *numneigh, **firstneigh;
double rcut, r_geom, delta, r_tri, r_particle, touch_distance, dt_crit;
int tri, particle;
Vector3d normal, x1, x2, x3, x4, x13, x23, x43, w, cp, x4cp, vnew, v_old;
;
Vector3d xi, x_center, dx;
Matrix2d C;
Vector2d w2d, rhs;
evdwl = 0.0;
ev_init(eflag, vflag);
tagint *mol = atom->molecule;
double **f = atom->f;
double **smd_data_9 = atom->smd_data_9;
double **x = atom->x;
double **x0 = atom->x0;
double **v = atom->v;
double *rmass = atom->rmass;
int *type = atom->type;
int nlocal = atom->nlocal;
double *radius = atom->contact_radius;
double rcutSq;
Vector3d offset;
int newton_pair = force->newton_pair;
int periodic = (domain->xperiodic || domain->yperiodic || domain->zperiodic);
inum = list->inum;
ilist = list->ilist;
numneigh = list->numneigh;
firstneigh = list->firstneigh;
int max_neighs = 0;
stable_time_increment = 1.0e22;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
itype = type[i];
jlist = firstneigh[i];
jnum = numneigh[i];
max_neighs = MAX(max_neighs, jnum);
for (jj = 0; jj < jnum; jj++) {
j = jlist[jj];
j &= NEIGHMASK;
jtype = type[j];
if ((mol[i] < 65535) && (mol[j] >= 65535)) {
particle = i;
tri = j;
} else if ((mol[j] < 65535) && (mol[i] >= 65535)) {
particle = j;
tri = i;
} else {
error->one(FLERR, "unknown case");
}
x_center(0) = x[tri][0];
x_center(1) = x[tri][1];
x_center(2) = x[tri][2];
x4(0) = x[particle][0];
x4(1) = x[particle][1];
x4(2) = x[particle][2];
dx = x_center - x4; if (periodic) {
domain->minimum_image(dx(0), dx(1), dx(2));
}
rsq = dx.squaredNorm();
r_tri = scale * radius[tri];
r_particle = scale * radius[particle];
rcut = r_tri + r_particle;
rcutSq = rcut * rcut;
if (rsq < rcutSq) {
normal(0) = x0[tri][0];
normal(1) = x0[tri][1];
normal(2) = x0[tri][2];
if (fabs(dx.dot(normal)) < radius[particle]) {
x1(0) = smd_data_9[tri][0];
x1(1) = smd_data_9[tri][1];
x1(2) = smd_data_9[tri][2];
x2(0) = smd_data_9[tri][3];
x2(1) = smd_data_9[tri][4];
x2(2) = smd_data_9[tri][5];
x3(0) = smd_data_9[tri][6];
x3(1) = smd_data_9[tri][7];
x3(2) = smd_data_9[tri][8];
PointTriangleDistance(x4, x1, x2, x3, cp, r);
x4cp = x4 - cp;
if (x4cp.dot(normal) < 0.0) {
normal *= -1.0;
}
if (r < 1.0 * radius[particle]) {
delta = radius[particle] - r; r_geom = radius[particle];
fpair = 1.066666667e0 * bulkmodulus[itype][jtype] * delta * sqrt(delta * r_geom);
dt_crit = 3.14 * sqrt(rmass[particle] / (fpair / delta));
stable_time_increment = MIN(stable_time_increment, dt_crit);
evdwl = r * fpair * 0.4e0 * delta;
fpair /= (r + 1.0e-2 * radius[particle]);
if (particle < nlocal) {
f[particle][0] += x4cp(0) * fpair;
f[particle][1] += x4cp(1) * fpair;
f[particle][2] += x4cp(2) * fpair;
}
if (tri < nlocal) {
f[tri][0] -= x4cp(0) * fpair;
f[tri][1] -= x4cp(1) * fpair;
f[tri][2] -= x4cp(2) * fpair;
}
if (evflag) {
ev_tally(i, j, nlocal, newton_pair, evdwl, 0.0, fpair, x4cp(0), x4cp(1), x4cp(2));
}
}
touch_distance = 1.0 * radius[particle];
if (r < touch_distance) {
normal = x4cp / r;
v_old(0) = v[particle][0];
v_old(1) = v[particle][1];
v_old(2) = v[particle][2];
if (v_old.dot(normal) < 0.0) {
vnew = 1.0 * (-2.0 * v_old.dot(normal) * normal + v_old);
v[particle][0] = vnew(0);
v[particle][1] = vnew(1);
v[particle][2] = vnew(2);
}
x[particle][0] = cp(0) + touch_distance * normal(0);
x[particle][1] = cp(1) + touch_distance * normal(1);
x[particle][2] = cp(2) + touch_distance * normal(2);
}
}
}
}
}
}
void PairTriSurf::allocate() {
allocated = 1;
int n = atom->ntypes;
memory->create(setflag, n + 1, n + 1, "pair:setflag");
for (int i = 1; i <= n; i++)
for (int j = i; j <= n; j++)
setflag[i][j] = 0;
memory->create(bulkmodulus, n + 1, n + 1, "pair:kspring");
memory->create(kn, n + 1, n + 1, "pair:kn");
memory->create(cutsq, n + 1, n + 1, "pair:cutsq");
onerad_dynamic = new double[n + 1];
onerad_frozen = new double[n + 1];
maxrad_dynamic = new double[n + 1];
maxrad_frozen = new double[n + 1];
}
void PairTriSurf::settings(int narg, char **arg) {
if (narg != 1)
error->all(FLERR, "Illegal number of args for pair_style smd/tri_surface");
scale = force->numeric(FLERR, arg[0]);
if (comm->me == 0) {
printf("\n>>========>>========>>========>>========>>========>>========>>========>>========\n");
printf("SMD/TRI_SURFACE CONTACT SETTINGS:\n");
printf("... effective contact radius is scaled by %f\n", scale);
printf(">>========>>========>>========>>========>>========>>========>>========>>========\n");
}
}
void PairTriSurf::coeff(int narg, char **arg) {
if (narg != 3)
error->all(FLERR, "Incorrect args for pair coefficients");
if (!allocated)
allocate();
int ilo, ihi, jlo, jhi;
force->bounds(FLERR,arg[0], atom->ntypes, ilo, ihi);
force->bounds(FLERR,arg[1], atom->ntypes, jlo, jhi);
double bulkmodulus_one = atof(arg[2]);
double kn_one = 0.0;
if (domain->dimension == 3) {
kn_one = (16. / 15.) * bulkmodulus_one; } else {
kn_one = 0.251856195 * (2. / 3.) * bulkmodulus_one; }
int count = 0;
for (int i = ilo; i <= ihi; i++) {
for (int j = MAX(jlo, i); j <= jhi; j++) {
bulkmodulus[i][j] = bulkmodulus_one;
kn[i][j] = kn_one;
setflag[i][j] = 1;
count++;
}
}
if (count == 0)
error->all(FLERR, "Incorrect args for pair coefficients");
}
double PairTriSurf::init_one(int i, int j) {
if (!allocated)
allocate();
if (setflag[i][j] == 0)
error->all(FLERR, "All pair coeffs are not set");
bulkmodulus[j][i] = bulkmodulus[i][j];
kn[j][i] = kn[i][j];
double cutoff = maxrad_dynamic[i] + maxrad_dynamic[j];
cutoff = MAX(cutoff, maxrad_frozen[i] + maxrad_dynamic[j]);
cutoff = MAX(cutoff, maxrad_dynamic[i] + maxrad_frozen[j]);
if (comm->me == 0) {
printf("cutoff for pair smd/smd/tri_surface = %f\n", cutoff);
}
return cutoff;
}
void PairTriSurf::init_style() {
int i;
if (!atom->contact_radius_flag)
error->all(FLERR, "Pair style smd/smd/tri_surface requires atom style with contact_radius");
int irequest = neighbor->request(this);
neighbor->requests[irequest]->size = 1;
for (i = 1; i <= atom->ntypes; i++)
onerad_dynamic[i] = onerad_frozen[i] = 0.0;
double *radius = atom->radius;
int *type = atom->type;
int nlocal = atom->nlocal;
for (i = 0; i < nlocal; i++) {
onerad_dynamic[type[i]] = MAX(onerad_dynamic[type[i]], radius[i]);
}
MPI_Allreduce(&onerad_dynamic[1], &maxrad_dynamic[1], atom->ntypes, MPI_DOUBLE, MPI_MAX, world);
MPI_Allreduce(&onerad_frozen[1], &maxrad_frozen[1], atom->ntypes, MPI_DOUBLE, MPI_MAX, world);
}
void PairTriSurf::init_list(int id, NeighList *ptr) {
if (id == 0)
list = ptr;
}
double PairTriSurf::memory_usage() {
return 0.0;
}
void PairTriSurf::PointTriangleDistance(const Vector3d sourcePosition, const Vector3d TRI0, const Vector3d TRI1,
const Vector3d TRI2, Vector3d &CP, double &dist) {
Vector3d edge0 = TRI1 - TRI0;
Vector3d edge1 = TRI2 - TRI0;
Vector3d v0 = TRI0 - sourcePosition;
double a = edge0.dot(edge0);
double b = edge0.dot(edge1);
double c = edge1.dot(edge1);
double d = edge0.dot(v0);
double e = edge1.dot(v0);
double det = a * c - b * b;
double s = b * e - c * d;
double t = b * d - a * e;
if (s + t < det) {
if (s < 0.f) {
if (t < 0.f) {
if (d < 0.f) {
s = clamp(-d / a, 0.f, 1.f);
t = 0.f;
} else {
s = 0.f;
t = clamp(-e / c, 0.f, 1.f);
}
} else {
s = 0.f;
t = clamp(-e / c, 0.f, 1.f);
}
} else if (t < 0.f) {
s = clamp(-d / a, 0.f, 1.f);
t = 0.f;
} else {
float invDet = 1.f / det;
s *= invDet;
t *= invDet;
}
} else {
if (s < 0.f) {
float tmp0 = b + d;
float tmp1 = c + e;
if (tmp1 > tmp0) {
float numer = tmp1 - tmp0;
float denom = a - 2 * b + c;
s = clamp(numer / denom, 0.f, 1.f);
t = 1 - s;
} else {
t = clamp(-e / c, 0.f, 1.f);
s = 0.f;
}
} else if (t < 0.f) {
if (a + d > b + e) {
float numer = c + e - b - d;
float denom = a - 2 * b + c;
s = clamp(numer / denom, 0.f, 1.f);
t = 1 - s;
} else {
s = clamp(-e / c, 0.f, 1.f);
t = 0.f;
}
} else {
float numer = c + e - b - d;
float denom = a - 2 * b + c;
s = clamp(numer / denom, 0.f, 1.f);
t = 1.f - s;
}
}
CP = TRI0 + s * edge0 + t * edge1;
dist = (CP - sourcePosition).norm();
}
double PairTriSurf::clamp(const double a, const double min, const double max) {
if (a < min) {
return min;
} else if (a > max) {
return max;
} else {
return a;
}
}
void *PairTriSurf::extract(const char *str, int &) {
if (strcmp(str, "smd/tri_surface/stable_time_increment_ptr") == 0) {
return (void *) &stable_time_increment;
}
return NULL;
}