mirror of https://github.com/lammps/lammps.git
git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@274 f3b2605a-c512-4ea7-a41b-209d697bcdaa
This commit is contained in:
parent
4b48dc819d
commit
8ed88c02de
|
@ -8,8 +8,6 @@ SHELL = /bin/sh
|
|||
CC = icc
|
||||
CCFLAGS = -O -I/home/sjplimp/tools/mpich/include \
|
||||
-I/home/sjplimp/tools/fftw/include -DFFT_FFTW -DGZIP
|
||||
#CCFLAGS = -O -restrict -I/home/sjplimp/tools/mpich/include \
|
||||
# -I/home/sjplimp/tools/fftw/include -DFFT_FFTW -DGZIP
|
||||
DEPFLAGS = -M
|
||||
LINK = icc
|
||||
LINKFLAGS = -O -L/home/sjplimp/tools/mpich/lib \
|
||||
|
|
|
@ -25,12 +25,7 @@ using namespace LAMMPS_NS;
|
|||
|
||||
/* ---------------------------------------------------------------------- */
|
||||
|
||||
PairEAMOpt::PairEAMOpt(LAMMPS *lmp) : PairEAM(lmp)
|
||||
{
|
||||
rhor_ = NULL;
|
||||
frho_ = NULL;
|
||||
z2r_ = NULL;
|
||||
}
|
||||
PairEAMOpt::PairEAMOpt(LAMMPS *lmp) : PairEAM(lmp) {}
|
||||
|
||||
/* ---------------------------------------------------------------------- */
|
||||
|
||||
|
|
|
@ -45,11 +45,7 @@ class PairEAMOpt : virtual public PairEAM {
|
|||
virtual ~PairEAMOpt() {}
|
||||
void compute(int, int);
|
||||
|
||||
protected:
|
||||
double* restrict rhor_;
|
||||
double* restrict frho_;
|
||||
double* restrict z2r_;
|
||||
|
||||
private:
|
||||
template < int EFLAG, int VFLAG, int NEWTON_PAIR > void eval();
|
||||
};
|
||||
|
||||
|
@ -75,8 +71,8 @@ void PairEAMOpt::eval()
|
|||
double _pad[3];
|
||||
} fast_gamma_t;
|
||||
|
||||
double** restrict f;
|
||||
double* restrict coeff;
|
||||
double** __restrict__ f;
|
||||
double* __restrict__ coeff;
|
||||
|
||||
// grow energy array if necessary
|
||||
|
||||
|
@ -94,28 +90,28 @@ void PairEAMOpt::eval()
|
|||
if (VFLAG == 2) f = update->f_pair;
|
||||
else f = atom->f;
|
||||
|
||||
double** restrict x = atom->x;
|
||||
int* restrict type = atom->type;
|
||||
double** __restrict__ x = atom->x;
|
||||
int* __restrict__ type = atom->type;
|
||||
int nlocal = atom->nlocal;
|
||||
|
||||
vec3_t* restrict xx = (vec3_t*)x[0];
|
||||
vec3_t* restrict ff = (vec3_t*)f[0];
|
||||
vec3_t* __restrict__ xx = (vec3_t*)x[0];
|
||||
vec3_t* __restrict__ ff = (vec3_t*)f[0];
|
||||
|
||||
double tmp_cutforcesq = cutforcesq;
|
||||
double tmp_rdr = rdr;
|
||||
int nr2 = nr-2;
|
||||
int nr1 = nr-1;
|
||||
|
||||
int** restrict firstneigh = neighbor->firstneigh;
|
||||
int* restrict num = neighbor->numneigh;
|
||||
int** __restrict__ firstneigh = neighbor->firstneigh;
|
||||
int* __restrict__ num = neighbor->numneigh;
|
||||
|
||||
int ntypes = atom->ntypes;
|
||||
int ntypes2 = ntypes*ntypes;
|
||||
|
||||
fast_alpha_t* restrict fast_alpha =
|
||||
fast_alpha_t* __restrict__ fast_alpha =
|
||||
(fast_alpha_t*) malloc(ntypes2*(nr+1)*sizeof(fast_alpha_t));
|
||||
for( int i = 0; i < ntypes; i++) for( int j = 0; j < ntypes; j++) {
|
||||
fast_alpha_t* restrict tab = &fast_alpha[i*ntypes*nr+j*nr];
|
||||
fast_alpha_t* __restrict__ tab = &fast_alpha[i*ntypes*nr+j*nr];
|
||||
for(int m = 1; m <= nr; m++) {
|
||||
tab[m].rhor0i = rhor_spline[type2rhor[i+1][j+1]][m][6];
|
||||
tab[m].rhor1i = rhor_spline[type2rhor[i+1][j+1]][m][5];
|
||||
|
@ -127,12 +123,12 @@ void PairEAMOpt::eval()
|
|||
tab[m].rhor3j = rhor_spline[type2rhor[j+1][i+1]][m][3];
|
||||
}
|
||||
}
|
||||
fast_alpha_t* restrict tabeight = fast_alpha;
|
||||
fast_alpha_t* __restrict__ tabeight = fast_alpha;
|
||||
|
||||
fast_gamma_t* restrict fast_gamma =
|
||||
fast_gamma_t* __restrict__ fast_gamma =
|
||||
(fast_gamma_t*) malloc(ntypes2*(nr+1)*sizeof(fast_gamma_t));
|
||||
for( int i = 0; i < ntypes; i++) for( int j = 0; j < ntypes; j++) {
|
||||
fast_gamma_t* restrict tab = &fast_gamma[i*ntypes*nr+j*nr];
|
||||
fast_gamma_t* __restrict__ tab = &fast_gamma[i*ntypes*nr+j*nr];
|
||||
for(int m = 1; m <= nr; m++) {
|
||||
tab[m].rhor4i = rhor_spline[type2rhor[i+1][j+1]][m][2];
|
||||
tab[m].rhor5i = rhor_spline[type2rhor[i+1][j+1]][m][1];
|
||||
|
@ -149,7 +145,7 @@ void PairEAMOpt::eval()
|
|||
tab[m].z2r6 = z2r_spline[type2z2r[i+1][j+1]][m][0];
|
||||
}
|
||||
}
|
||||
fast_gamma_t* restrict tabss = fast_gamma;
|
||||
fast_gamma_t* __restrict__ tabss = fast_gamma;
|
||||
|
||||
// zero out density
|
||||
|
||||
|
@ -166,12 +162,12 @@ void PairEAMOpt::eval()
|
|||
double ytmp = xx[i].y;
|
||||
double ztmp = xx[i].z;
|
||||
int itype = type[i] - 1;
|
||||
int* restrict neighs = firstneigh[i];
|
||||
int* __restrict__ neighs = firstneigh[i];
|
||||
int numneigh = num[i];
|
||||
|
||||
double tmprho = rho[i];
|
||||
|
||||
fast_alpha_t* restrict tabeighti = &tabeight[itype*ntypes*nr];
|
||||
fast_alpha_t* __restrict__ tabeighti = &tabeight[itype*ntypes*nr];
|
||||
for (int k = 0; k < numneigh; k++) {
|
||||
int j = neighs[k];
|
||||
|
||||
|
@ -238,7 +234,7 @@ void PairEAMOpt::eval()
|
|||
double ztmp = xx[i].z;
|
||||
int itype = type[i];
|
||||
int itype1 = type[i] - 1;
|
||||
int* restrict neighs = firstneigh[i];
|
||||
int* __restrict__ neighs = firstneigh[i];
|
||||
int numneigh = num[i];
|
||||
|
||||
double tmpfx = 0.0;
|
||||
|
@ -247,7 +243,7 @@ void PairEAMOpt::eval()
|
|||
|
||||
double fpi = fp[i];
|
||||
|
||||
fast_gamma_t* restrict tabssi = &tabss[itype1*ntypes*nr];
|
||||
fast_gamma_t* __restrict__ tabssi = &tabss[itype1*ntypes*nr];
|
||||
for (int k = 0; k < numneigh; k++) {
|
||||
int j = neighs[k];
|
||||
|
||||
|
|
|
@ -75,8 +75,8 @@ void PairLJCharmmCoulLongOpt::eval()
|
|||
double grij,expm2,prefactor,t,erfc;
|
||||
double factor,phicoul,philj,switch1,switch2;
|
||||
|
||||
int* restrict neighs;
|
||||
double** restrict f;
|
||||
int* __restrict__ neighs;
|
||||
double** __restrict__ f;
|
||||
|
||||
float rsq;
|
||||
int *int_rsq = (int *) &rsq;
|
||||
|
@ -87,19 +87,19 @@ void PairLJCharmmCoulLongOpt::eval()
|
|||
if (VFLAG == 2) f = update->f_pair;
|
||||
else f = atom->f;
|
||||
|
||||
double** restrict x = atom->x;
|
||||
double* restrict q = atom->q;
|
||||
int* restrict type = atom->type;
|
||||
double** __restrict__ x = atom->x;
|
||||
double* __restrict__ q = atom->q;
|
||||
int* __restrict__ type = atom->type;
|
||||
int nlocal = atom->nlocal;
|
||||
int nall = atom->nlocal + atom->nghost;
|
||||
double* restrict special_coul = force->special_coul;
|
||||
double* restrict special_lj = force->special_lj;
|
||||
double* __restrict__ special_coul = force->special_coul;
|
||||
double* __restrict__ special_lj = force->special_lj;
|
||||
double qqrd2e = force->qqrd2e;
|
||||
int** restrict firstneigh = neighbor->firstneigh;
|
||||
int* restrict num = neighbor->numneigh;
|
||||
int** __restrict__ firstneigh = neighbor->firstneigh;
|
||||
int* __restrict__ num = neighbor->numneigh;
|
||||
|
||||
vec3_t* restrict xx = (vec3_t*)x[0];
|
||||
vec3_t* restrict ff = (vec3_t*)f[0];
|
||||
vec3_t* __restrict__ xx = (vec3_t*)x[0];
|
||||
vec3_t* __restrict__ ff = (vec3_t*)f[0];
|
||||
|
||||
int ntypes = atom->ntypes;
|
||||
int ntypes2 = ntypes*ntypes;
|
||||
|
@ -107,7 +107,7 @@ void PairLJCharmmCoulLongOpt::eval()
|
|||
double tmp_coef1 = 1.0/denom_lj;
|
||||
double tmp_coef2 = cut_ljsq - 3.0*cut_lj_innersq;
|
||||
|
||||
fast_alpha_t* restrict fast_alpha =
|
||||
fast_alpha_t* __restrict__ fast_alpha =
|
||||
(fast_alpha_t*)malloc(ntypes2*sizeof(fast_alpha_t));
|
||||
for( int i = 0; i < ntypes; i++) for( int j = 0; j < ntypes; j++) {
|
||||
fast_alpha_t& a = fast_alpha[i*ntypes+j];
|
||||
|
@ -117,7 +117,7 @@ void PairLJCharmmCoulLongOpt::eval()
|
|||
a.lj3 = lj3[i+1][j+1];
|
||||
a.lj4 = lj4[i+1][j+1];
|
||||
}
|
||||
fast_alpha_t* restrict tabsix = fast_alpha;
|
||||
fast_alpha_t* __restrict__ tabsix = fast_alpha;
|
||||
|
||||
// loop over neighbors of my atoms
|
||||
|
||||
|
@ -127,14 +127,14 @@ void PairLJCharmmCoulLongOpt::eval()
|
|||
double ytmp = xx[i].y;
|
||||
double ztmp = xx[i].z;
|
||||
itype = type[i] - 1;
|
||||
int* restrict neighs = firstneigh[i];
|
||||
int* __restrict__ neighs = firstneigh[i];
|
||||
int numneigh = num[i];
|
||||
|
||||
double tmpfx = 0.0;
|
||||
double tmpfy = 0.0;
|
||||
double tmpfz = 0.0;
|
||||
|
||||
fast_alpha_t* restrict tabsixi = (fast_alpha_t*) &tabsix[itype*ntypes];
|
||||
fast_alpha_t* __restrict__ tabsixi = (fast_alpha_t*) &tabsix[itype*ntypes];
|
||||
for (k = 0; k < numneigh; k++) {
|
||||
j = neighs[k];
|
||||
|
||||
|
|
|
@ -49,8 +49,8 @@ void PairLJCutOpt::eval()
|
|||
double _pad[2];
|
||||
} fast_alpha_t;
|
||||
|
||||
int* restrict neighs;
|
||||
double** restrict f;
|
||||
int* __restrict__ neighs;
|
||||
double** __restrict__ f;
|
||||
|
||||
eng_vdwl = 0.0;
|
||||
if (VFLAG) for (int i = 0; i < 6; i++) virial[i] = 0.0;
|
||||
|
@ -58,21 +58,21 @@ void PairLJCutOpt::eval()
|
|||
if (VFLAG == 2) f = update->f_pair;
|
||||
else f = atom->f;
|
||||
|
||||
double** restrict x = atom->x;
|
||||
int* restrict type = atom->type;
|
||||
double** __restrict__ x = atom->x;
|
||||
int* __restrict__ type = atom->type;
|
||||
int nlocal = atom->nlocal;
|
||||
int nall = atom->nlocal + atom->nghost;
|
||||
double* restrict special_lj = force->special_lj;
|
||||
int** restrict firstneigh = neighbor->firstneigh;
|
||||
int* restrict num = neighbor->numneigh;
|
||||
double* __restrict__ special_lj = force->special_lj;
|
||||
int** __restrict__ firstneigh = neighbor->firstneigh;
|
||||
int* __restrict__ num = neighbor->numneigh;
|
||||
|
||||
vec3_t* restrict xx = (vec3_t*)x[0];
|
||||
vec3_t* restrict ff = (vec3_t*)f[0];
|
||||
vec3_t* __restrict__ xx = (vec3_t*)x[0];
|
||||
vec3_t* __restrict__ ff = (vec3_t*)f[0];
|
||||
|
||||
int ntypes = atom->ntypes;
|
||||
int ntypes2 = ntypes*ntypes;
|
||||
|
||||
fast_alpha_t* restrict fast_alpha =
|
||||
fast_alpha_t* __restrict__ fast_alpha =
|
||||
(fast_alpha_t*) malloc(ntypes2*sizeof(fast_alpha_t));
|
||||
for( int i = 0; i < ntypes; i++) for( int j = 0; j < ntypes; j++) {
|
||||
fast_alpha_t& a = fast_alpha[i*ntypes+j];
|
||||
|
@ -83,7 +83,7 @@ void PairLJCutOpt::eval()
|
|||
a.lj4 = lj4[i+1][j+1];
|
||||
a.offset = offset[i+1][j+1];
|
||||
}
|
||||
fast_alpha_t* restrict tabsix = fast_alpha;
|
||||
fast_alpha_t* __restrict__ tabsix = fast_alpha;
|
||||
|
||||
// loop over neighbors of my atoms
|
||||
|
||||
|
@ -92,14 +92,14 @@ void PairLJCutOpt::eval()
|
|||
double ytmp = xx[i].y;
|
||||
double ztmp = xx[i].z;
|
||||
int itype = type[i] - 1;
|
||||
int* restrict neighs = firstneigh[i];
|
||||
int* __restrict__ neighs = firstneigh[i];
|
||||
int numneigh = num[i];
|
||||
|
||||
double tmpfx = 0.0;
|
||||
double tmpfy = 0.0;
|
||||
double tmpfz = 0.0;
|
||||
|
||||
fast_alpha_t* restrict tabsixi = (fast_alpha_t*)&tabsix[itype*ntypes];
|
||||
fast_alpha_t* __restrict__ tabsixi = (fast_alpha_t*)&tabsix[itype*ntypes];
|
||||
|
||||
for (int k = 0; k < numneigh; k++) {
|
||||
int j = neighs[k];
|
||||
|
|
|
@ -54,8 +54,8 @@ void PairMorseOpt::eval()
|
|||
double _pad[2];
|
||||
} fast_alpha_t;
|
||||
|
||||
int* restrict neighs;
|
||||
double** restrict f;
|
||||
int* __restrict__ neighs;
|
||||
double** __restrict__ f;
|
||||
|
||||
eng_vdwl = 0.0;
|
||||
if (VFLAG) for (int i = 0; i < 6; i++) virial[i] = 0.0;
|
||||
|
@ -63,22 +63,22 @@ void PairMorseOpt::eval()
|
|||
if (VFLAG == 2) f = update->f_pair;
|
||||
else f = atom->f;
|
||||
|
||||
double** restrict x = atom->x;
|
||||
int* restrict type = atom->type;
|
||||
double** __restrict__ x = atom->x;
|
||||
int* __restrict__ type = atom->type;
|
||||
int nlocal = atom->nlocal;
|
||||
int nall = atom->nlocal + atom->nghost;
|
||||
double* restrict special_lj = force->special_lj;
|
||||
double* __restrict__ special_lj = force->special_lj;
|
||||
|
||||
int** restrict firstneigh = neighbor->firstneigh;
|
||||
int* restrict num = neighbor->numneigh;
|
||||
int** __restrict__ firstneigh = neighbor->firstneigh;
|
||||
int* __restrict__ num = neighbor->numneigh;
|
||||
|
||||
vec3_t* restrict xx = (vec3_t*)x[0];
|
||||
vec3_t* restrict ff = (vec3_t*)f[0];
|
||||
vec3_t* __restrict__ xx = (vec3_t*)x[0];
|
||||
vec3_t* __restrict__ ff = (vec3_t*)f[0];
|
||||
|
||||
int ntypes = atom->ntypes;
|
||||
int ntypes2 = ntypes*ntypes;
|
||||
|
||||
fast_alpha_t* restrict fast_alpha =
|
||||
fast_alpha_t* __restrict__ fast_alpha =
|
||||
(fast_alpha_t*) malloc(ntypes2*sizeof(fast_alpha_t));
|
||||
for( int i = 0; i < ntypes; i++) for( int j = 0; j < ntypes; j++) {
|
||||
fast_alpha_t& a = fast_alpha[i*ntypes+j];
|
||||
|
@ -89,7 +89,7 @@ void PairMorseOpt::eval()
|
|||
a.d0 = d0[i+1][j+1];
|
||||
a.offset = offset[i+1][j+1];
|
||||
}
|
||||
fast_alpha_t* restrict tabsix = fast_alpha;
|
||||
fast_alpha_t* __restrict__ tabsix = fast_alpha;
|
||||
|
||||
// loop over neighbors of my atoms
|
||||
|
||||
|
@ -98,14 +98,14 @@ void PairMorseOpt::eval()
|
|||
double ytmp = xx[i].y;
|
||||
double ztmp = xx[i].z;
|
||||
int itype = type[i] - 1;
|
||||
int* restrict neighs = firstneigh[i];
|
||||
int* __restrict__ neighs = firstneigh[i];
|
||||
int numneigh = num[i];
|
||||
|
||||
double tmpfx = 0.0;
|
||||
double tmpfy = 0.0;
|
||||
double tmpfz = 0.0;
|
||||
|
||||
fast_alpha_t* restrict tabsixi = (fast_alpha_t*)&tabsix[itype*ntypes];
|
||||
fast_alpha_t* __restrict__ tabsixi = (fast_alpha_t*)&tabsix[itype*ntypes];
|
||||
|
||||
for (int k = 0; k < numneigh; k++) {
|
||||
int j = neighs[k];
|
||||
|
|
|
@ -0,0 +1,30 @@
|
|||
/* ----------------------------------------------------------------------
|
||||
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
||||
http://lammps.sandia.gov, Sandia National Laboratories
|
||||
Steve Plimpton, sjplimp@sandia.gov
|
||||
|
||||
Copyright (2003) Sandia Corporation. Under the terms of Contract
|
||||
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
|
||||
certain rights in this software. This software is distributed under
|
||||
the GNU General Public License.
|
||||
|
||||
See the README file in the top-level LAMMPS directory.
|
||||
------------------------------------------------------------------------- */
|
||||
|
||||
#ifdef PairInclude
|
||||
#include "pair_eam_opt.h"
|
||||
#include "pair_eam_alloy_opt.h"
|
||||
#include "pair_eam_fs_opt.h"
|
||||
#include "pair_lj_charmm_coul_long_opt.h"
|
||||
#include "pair_lj_cut_opt.h"
|
||||
#include "pair_morse_opt.h"
|
||||
#endif
|
||||
|
||||
#ifdef PairClass
|
||||
PairStyle(eam/opt,PairEAMOpt)
|
||||
PairStyle(eam/alloy/opt,PairEAMAlloyOpt)
|
||||
PairStyle(eam/fs/opt,PairEAMFSOpt)
|
||||
PairStyle(lj/cut/opt,PairLJCutOpt)
|
||||
PairStyle(lj/charmm/coul/long/opt,PairLJCharmmCoulLongOpt)
|
||||
PairStyle(morse/opt,PairMorseOpt)
|
||||
#endif
|
Loading…
Reference in New Issue