1 // clang-format off
2 /* ----------------------------------------------------------------------
3 LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
4 https://www.lammps.org/, Sandia National Laboratories
5 Steve Plimpton, sjplimp@sandia.gov
6
7 Copyright (2003) Sandia Corporation. Under the terms of Contract
8 DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
9 certain rights in this software. This software is distributed under
10 the GNU General Public License.
11
12 See the README file in the top-level LAMMPS directory.
13 ------------------------------------------------------------------------- */
14
15 /* ----------------------------------------------------------------------
16 Contributing author: Jan Los
17 ------------------------------------------------------------------------- */
18
19 #include "pair_extep.h"
20
21 #include "atom.h"
22 #include "comm.h"
23 #include "error.h"
24 #include "force.h"
25 #include "math_const.h"
26 #include "math_extra.h"
27 #include "memory.h"
28 #include "my_page.h"
29 #include "neigh_list.h"
30 #include "neigh_request.h"
31 #include "neighbor.h"
32
33 #include <cmath>
34 #include <cstring>
35 #include <cctype>
36
37 using namespace LAMMPS_NS;
38 using namespace MathConst;
39 using namespace MathExtra;
40
41 #define MAXLINE 1024
42 #define DELTA 4
43 #define PGDELTA 1
44
45 /* ---------------------------------------------------------------------- */
46
PairExTeP(LAMMPS * lmp)47 PairExTeP::PairExTeP(LAMMPS *lmp) : Pair(lmp)
48 {
49 single_enable = 0;
50 restartinfo = 0;
51 one_coeff = 1;
52 manybody_flag = 1;
53 centroidstressflag = CENTROID_NOTAVAIL;
54 ghostneigh = 1;
55
56 params = nullptr;
57
58 maxlocal = 0;
59 SR_numneigh = nullptr;
60 SR_firstneigh = nullptr;
61 ipage = nullptr;
62 pgsize = oneatom = 0;
63
64 Nt = nullptr;
65 Nd = nullptr;
66 }
67
68 /* ----------------------------------------------------------------------
69 check if allocated, since class can be destructed when incomplete
70 ------------------------------------------------------------------------- */
71
~PairExTeP()72 PairExTeP::~PairExTeP()
73 {
74 memory->destroy(params);
75 memory->destroy(elem3param);
76
77 memory->destroy(SR_numneigh);
78 memory->sfree(SR_firstneigh);
79 delete [] ipage;
80 memory->destroy(Nt);
81 memory->destroy(Nd);
82
83 if (allocated) {
84 memory->destroy(setflag);
85 memory->destroy(cutsq);
86 memory->destroy(cutghost);
87 }
88 }
89
90 /* ----------------------------------------------------------------------
91 create SR neighbor list from main neighbor list
92 SR neighbor list stores neighbors of ghost atoms
93 ------------------------------------------------------------------------- */
94
SR_neigh()95 void PairExTeP::SR_neigh()
96 {
97 int i,j,ii,jj,n,allnum,jnum,itype,jtype,iparam_ij;
98 double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
99 int *ilist,*jlist,*numneigh,**firstneigh;
100 int *neighptr;
101
102 double **x = atom->x;
103 int *type = atom->type;
104
105 if (atom->nmax > maxlocal) { // ensure there is enough space
106 maxlocal = atom->nmax; // for atoms and ghosts allocated
107 memory->destroy(SR_numneigh);
108 memory->sfree(SR_firstneigh);
109 memory->destroy(Nt);
110 memory->destroy(Nd);
111 memory->create(SR_numneigh,maxlocal,"ExTeP:numneigh");
112 SR_firstneigh = (int **) memory->smalloc(maxlocal*sizeof(int *),
113 "ExTeP:firstneigh");
114 memory->create(Nt,maxlocal,"ExTeP:Nt");
115 memory->create(Nd,maxlocal,"ExTeP:Nd");
116 }
117
118 allnum = list->inum + list->gnum;
119 ilist = list->ilist;
120 numneigh = list->numneigh;
121 firstneigh = list->firstneigh;
122
123 // store all SR neighs of owned and ghost atoms
124 // scan full neighbor list of I
125
126 ipage->reset();
127
128 for (ii = 0; ii < allnum; ii++) {
129 i = ilist[ii];
130 itype=map[type[i]];
131
132 n = 0;
133 neighptr = ipage->vget();
134
135 xtmp = x[i][0];
136 ytmp = x[i][1];
137 ztmp = x[i][2];
138
139 Nt[i] = 0.0;
140 Nd[i] = 0.0;
141
142 jlist = firstneigh[i];
143 jnum = numneigh[i];
144
145 for (jj = 0; jj < jnum; jj++) {
146 j = jlist[jj];
147 j &= NEIGHMASK;
148 delx = xtmp - x[j][0];
149 dely = ytmp - x[j][1];
150 delz = ztmp - x[j][2];
151 rsq = delx*delx + dely*dely + delz*delz;
152
153 jtype=map[type[j]];
154 iparam_ij = elem3param[itype][jtype][jtype];
155
156 if (rsq < params[iparam_ij].cutsq) {
157 neighptr[n++] = j;
158 double tmp_fc = ters_fc(sqrt(rsq),¶ms[iparam_ij]);
159 Nt[i] += tmp_fc;
160 if (itype!=jtype) {
161 Nd[i] += tmp_fc;
162 }
163 }
164 }
165 //printf("SR_neigh : N[%d] = %f\n",i,N[i]);
166
167 ipage->vgot(n);
168 if (ipage->status())
169 error->one(FLERR,"Neighbor list overflow, boost neigh_modify one");
170 }
171 }
172
173
174 /* ---------------------------------------------------------------------- */
175
compute(int eflag,int vflag)176 void PairExTeP::compute(int eflag, int vflag)
177 {
178 int i,j,k,ii,jj,kk,inum,jnum;
179 int itype,jtype,ktype,iparam_ij,iparam_ijk;
180 tagint itag,jtag;
181 double xtmp,ytmp,ztmp,delx,dely,delz,evdwl,fpair;
182 double rsq,rsq1,rsq2,r2;
183 double delr1[3],delr2[3],fi[3],fj[3],fk[3];
184 double zeta_ij,prefactor;
185 int *ilist,*jlist,*numneigh,**firstneigh;
186
187 evdwl = 0.0;
188 ev_init(eflag,vflag);
189
190 SR_neigh();
191
192 double **x = atom->x;
193 double **f = atom->f;
194 tagint *tag = atom->tag;
195 int *type = atom->type;
196 int nlocal = atom->nlocal;
197 int newton_pair = force->newton_pair;
198
199 inum = list->inum;
200 ilist = list->ilist;
201 numneigh = list->numneigh;
202 firstneigh = list->firstneigh;
203
204 // loop over full neighbor list of my atoms
205
206 for (ii = 0; ii < inum; ii++) {
207 i = ilist[ii];
208 itag = tag[i];
209 itype = map[type[i]];
210 xtmp = x[i][0];
211 ytmp = x[i][1];
212 ztmp = x[i][2];
213
214 // two-body interactions, skip half of them
215
216 jlist = firstneigh[i];
217 jnum = numneigh[i];
218
219 for (jj = 0; jj < jnum; jj++) {
220 j = jlist[jj];
221 j &= NEIGHMASK;
222 jtag = tag[j];
223
224 if (itag > jtag) {
225 if ((itag+jtag) % 2 == 0) continue;
226 } else if (itag < jtag) {
227 if ((itag+jtag) % 2 == 1) continue;
228 } else {
229 if (x[j][2] < x[i][2]) continue;
230 if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
231 if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp) continue;
232 }
233
234 jtype = map[type[j]];
235
236 delx = xtmp - x[j][0];
237 dely = ytmp - x[j][1];
238 delz = ztmp - x[j][2];
239 rsq = delx*delx + dely*dely + delz*delz;
240
241 iparam_ij = elem3param[itype][jtype][jtype];
242 if (rsq > params[iparam_ij].cutsq) continue;
243
244 repulsive(¶ms[iparam_ij],rsq,fpair,eflag,evdwl);
245
246 f[i][0] += delx*fpair;
247 f[i][1] += dely*fpair;
248 f[i][2] += delz*fpair;
249 f[j][0] -= delx*fpair;
250 f[j][1] -= dely*fpair;
251 f[j][2] -= delz*fpair;
252
253 if (evflag) ev_tally(i,j,nlocal,newton_pair,
254 evdwl,0.0,fpair,delx,dely,delz);
255 }
256
257 // three-body interactions -(bij + Fcorrection) * fA
258 // skip immediately if I-J is not within cutoff
259
260 for (jj = 0; jj < jnum; jj++) {
261 j = jlist[jj];
262 j &= NEIGHMASK;
263 jtag = tag[j];
264 jtype = map[type[j]];
265 iparam_ij = elem3param[itype][jtype][jtype];
266
267 delr1[0] = x[j][0] - xtmp;
268 delr1[1] = x[j][1] - ytmp;
269 delr1[2] = x[j][2] - ztmp;
270 rsq1 = delr1[0]*delr1[0] + delr1[1]*delr1[1] + delr1[2]*delr1[2];
271 if (rsq1 > params[iparam_ij].cutsq) continue;
272
273 // accumulate bondorder zeta for each i-j interaction via loop over k
274
275 zeta_ij = 0.0;
276
277 /* F_IJ (1) */
278 // compute correction to energy and forces
279 // dE/dr = -Fij(Zi,Zj) dV/dr
280 // - dFij/dZi dZi/dr V
281 // (conjugate term is computed when j is a central atom)
282
283 double FXY, dFXY_dNdij, dFXY_dNdji, fa, fa_d, deng, fpair;
284 double Ntij = Nt[i];
285 double Ndij = Nd[i];
286 double Ntji = Nt[j];
287 double Ndji = Nd[j];
288 double r = sqrt(rsq1);
289 double fc_ij = ters_fc(r,¶ms[iparam_ij]);
290
291 Ntij -= fc_ij;
292 Ntji -= fc_ij;
293 if (jtype!=itype) {
294 Ndij -= fc_ij;
295 Ndji -= fc_ij;
296 }
297 if (Ntij<0) { Ntij=0.; }
298 if (Ndij<0) { Ndij=0.; }
299 if (Ntji<0) { Ntji=0.; }
300 if (Ndji<0) { Ndji=0.; }
301 FXY = F_corr(itype, jtype, Ndij, Ndji, &dFXY_dNdij, &dFXY_dNdji);
302
303 // envelop functions
304 double fenv, dfenv_ij;
305 fenv = envelop_function(Ntij, Ntji, &dfenv_ij);
306 //
307 double Fc = fenv * FXY;
308 double dFc_dNtij = dfenv_ij * FXY;
309 double dFc_dNdij = fenv * dFXY_dNdij;
310
311 fa = ters_fa(r,¶ms[iparam_ij]);
312 fa_d = ters_fa_d(r,¶ms[iparam_ij]);
313 deng = 0.5 * fa * Fc;
314 fpair = 0.5 * fa_d * Fc / r;
315
316 f[i][0] += delr1[0]*fpair;
317 f[i][1] += delr1[1]*fpair;
318 f[i][2] += delr1[2]*fpair;
319 f[j][0] -= delr1[0]*fpair;
320 f[j][1] -= delr1[1]*fpair;
321 f[j][2] -= delr1[2]*fpair;
322
323 if (evflag) ev_tally(i,j,nlocal,newton_pair,
324 deng,0.0,-fpair,-delr1[0],-delr1[1],-delr1[2]);
325 /* END F_IJ (1) */
326
327 for (kk = 0; kk < jnum; kk++) {
328 if (jj == kk) continue;
329 k = jlist[kk];
330 k &= NEIGHMASK;
331 ktype = map[type[k]];
332 iparam_ijk = elem3param[itype][jtype][ktype];
333
334 delr2[0] = x[k][0] - xtmp;
335 delr2[1] = x[k][1] - ytmp;
336 delr2[2] = x[k][2] - ztmp;
337 rsq2 = delr2[0]*delr2[0] + delr2[1]*delr2[1] + delr2[2]*delr2[2];
338 if (rsq2 > params[iparam_ijk].cutsq) continue;
339
340 r2 = sqrt(rsq2);
341
342 zeta_ij += zeta(¶ms[iparam_ijk],r,r2,delr1,delr2);
343
344 /* F_IJ (2) */
345 // compute force components due to spline derivatives
346 // uses only the part with FXY_x (FXY_y is done when i and j are inversed)
347 int iparam_ik = elem3param[itype][ktype][0];
348 double fc_ik_d = ters_fc_d(r2,¶ms[iparam_ik]);
349 double fc_prefac_ik_0 = 1.0 * fc_ik_d * fa / r2;
350 double fc_prefac_ik = dFc_dNtij * fc_prefac_ik_0;
351 f[i][0] += fc_prefac_ik * delr2[0];
352 f[i][1] += fc_prefac_ik * delr2[1];
353 f[i][2] += fc_prefac_ik * delr2[2];
354 f[k][0] -= fc_prefac_ik * delr2[0];
355 f[k][1] -= fc_prefac_ik * delr2[1];
356 f[k][2] -= fc_prefac_ik * delr2[2];
357 if (vflag_either) v_tally2(i,k,-fc_prefac_ik,delr2);
358 if (itype != ktype) {
359 fc_prefac_ik = dFc_dNdij * fc_prefac_ik_0;
360 f[i][0] += fc_prefac_ik * delr2[0];
361 f[i][1] += fc_prefac_ik * delr2[1];
362 f[i][2] += fc_prefac_ik * delr2[2];
363 f[k][0] -= fc_prefac_ik * delr2[0];
364 f[k][1] -= fc_prefac_ik * delr2[1];
365 f[k][2] -= fc_prefac_ik * delr2[2];
366 if (vflag_either) v_tally2(i,k,-fc_prefac_ik,delr2);
367 }
368 /* END F_IJ (2) */
369
370 }
371
372 // pairwise force due to zeta
373
374 force_zeta(¶ms[iparam_ij],r,zeta_ij,fpair,prefactor,eflag,evdwl);
375
376 f[i][0] += delr1[0]*fpair;
377 f[i][1] += delr1[1]*fpair;
378 f[i][2] += delr1[2]*fpair;
379 f[j][0] -= delr1[0]*fpair;
380 f[j][1] -= delr1[1]*fpair;
381 f[j][2] -= delr1[2]*fpair;
382
383 if (evflag) ev_tally(i,j,nlocal,newton_pair,
384 evdwl,0.0,-fpair,-delr1[0],-delr1[1],-delr1[2]);
385
386 // attractive term via loop over k
387
388 for (kk = 0; kk < jnum; kk++) {
389 if (jj == kk) continue;
390 k = jlist[kk];
391 k &= NEIGHMASK;
392 ktype = map[type[k]];
393 iparam_ijk = elem3param[itype][jtype][ktype];
394
395 delr2[0] = x[k][0] - xtmp;
396 delr2[1] = x[k][1] - ytmp;
397 delr2[2] = x[k][2] - ztmp;
398 rsq2 = delr2[0]*delr2[0] + delr2[1]*delr2[1] + delr2[2]*delr2[2];
399 if (rsq2 > params[iparam_ijk].cutsq) continue;
400
401 attractive(¶ms[iparam_ijk],prefactor,
402 rsq1,rsq2,delr1,delr2,fi,fj,fk);
403
404
405 f[i][0] += fi[0];
406 f[i][1] += fi[1];
407 f[i][2] += fi[2];
408 f[j][0] += fj[0];
409 f[j][1] += fj[1];
410 f[j][2] += fj[2];
411 f[k][0] += fk[0];
412 f[k][1] += fk[1];
413 f[k][2] += fk[2];
414
415 if (vflag_either) v_tally3(i,j,k,fj,fk,delr1,delr2);
416 }
417 }
418 }
419
420 if (vflag_fdotr) virial_fdotr_compute();
421 }
422
423 /* ---------------------------------------------------------------------- */
424
allocate()425 void PairExTeP::allocate()
426 {
427 allocated = 1;
428 int n = atom->ntypes;
429
430 memory->create(setflag,n+1,n+1,"pair:setflag");
431 memory->create(cutsq,n+1,n+1,"pair:cutsq");
432 memory->create(cutghost,n+1,n+1,"pair:cutghost");
433
434 map = new int[n+1];
435 }
436
437 /* ----------------------------------------------------------------------
438 global settings
439 ------------------------------------------------------------------------- */
440
settings(int narg,char **)441 void PairExTeP::settings(int narg, char **/*arg*/)
442 {
443 if (narg != 0) error->all(FLERR,"Illegal pair_style command");
444 }
445
446 /* ----------------------------------------------------------------------
447 set coeffs for one or more type pairs
448 ------------------------------------------------------------------------- */
449
coeff(int narg,char ** arg)450 void PairExTeP::coeff(int narg, char **arg)
451 {
452 if (!allocated) allocate();
453
454 map_element2type(narg-3,arg+3);
455
456 // read potential file and initialize potential parameters
457
458 read_file(arg[2]);
459 spline_init();
460 setup();
461 }
462
463 /* ----------------------------------------------------------------------
464 init specific to this pair style
465 ------------------------------------------------------------------------- */
466
init_style()467 void PairExTeP::init_style()
468 {
469 if (atom->tag_enable == 0)
470 error->all(FLERR,"Pair style ExTeP requires atom IDs");
471 if (force->newton_pair == 0)
472 error->all(FLERR,"Pair style ExTeP requires newton pair on");
473
474 // need a full neighbor list
475
476 int irequest = neighbor->request(this);
477 neighbor->requests[irequest]->half = 0;
478 neighbor->requests[irequest]->full = 1;
479
480 // including neighbors of ghosts
481 neighbor->requests[irequest]->ghost = 1;
482
483 // create pages if first time or if neighbor pgsize/oneatom has changed
484
485 int create = 0;
486 if (ipage == nullptr) create = 1;
487 if (pgsize != neighbor->pgsize) create = 1;
488 if (oneatom != neighbor->oneatom) create = 1;
489
490 if (create) {
491 delete [] ipage;
492 pgsize = neighbor->pgsize;
493 oneatom = neighbor->oneatom;
494
495 int nmypage= comm->nthreads;
496 ipage = new MyPage<int>[nmypage];
497 for (int i = 0; i < nmypage; i++)
498 ipage[i].init(oneatom,pgsize,PGDELTA);
499 }
500 }
501
502 /* ----------------------------------------------------------------------
503 init for one type pair i,j and corresponding j,i
504 ------------------------------------------------------------------------- */
505
init_one(int i,int j)506 double PairExTeP::init_one(int i, int j)
507 {
508 if (setflag[i][j] == 0) error->all(FLERR,"All pair coeffs are not set");
509
510 cutghost[i][j] = cutmax ;
511 cutghost[j][i] = cutghost[i][j];
512
513 return cutmax;
514 }
515
516 /* ---------------------------------------------------------------------- */
517
read_file(char * file)518 void PairExTeP::read_file(char *file)
519 {
520 int params_per_line = 17;
521 char **words = new char*[params_per_line+1];
522
523 memory->sfree(params);
524 params = nullptr;
525 nparams = maxparam = 0;
526
527 // open file on proc 0
528
529 FILE *fp;
530 if (comm->me == 0) {
531 fp = utils::open_potential(file,lmp,nullptr);
532 if (fp == nullptr)
533 error->one(FLERR,"Cannot open ExTeP potential file {}: {}",file,utils::getsyserror());
534 }
535
536 // read each line out of file, skipping blank lines or leading '#'
537 // store line of params if all 3 element tags are in element list
538
539 int n,nwords,ielement,jelement,kelement;
540 char line[MAXLINE],*ptr;
541 int eof = 0;
542
543 while (1) {
544 if (comm->me == 0) {
545 ptr = fgets(line,MAXLINE,fp);
546 if (ptr == nullptr) {
547 eof = 1;
548 fclose(fp);
549 } else n = strlen(line) + 1;
550 }
551 MPI_Bcast(&eof,1,MPI_INT,0,world);
552 if (eof) break;
553 MPI_Bcast(&n,1,MPI_INT,0,world);
554 MPI_Bcast(line,n,MPI_CHAR,0,world);
555
556 // strip comment, skip line if blank
557
558 if ((ptr = strchr(line,'#'))) *ptr = '\0';
559 nwords = utils::count_words(line);
560 if (nwords == 0) continue;
561
562 // concatenate additional lines until have params_per_line words
563
564 while (nwords < params_per_line) {
565 n = strlen(line);
566 if (comm->me == 0) {
567 ptr = fgets(&line[n],MAXLINE-n,fp);
568 if (ptr == nullptr) {
569 eof = 1;
570 fclose(fp);
571 } else n = strlen(line) + 1;
572 }
573 MPI_Bcast(&eof,1,MPI_INT,0,world);
574 if (eof) break;
575 MPI_Bcast(&n,1,MPI_INT,0,world);
576 MPI_Bcast(line,n,MPI_CHAR,0,world);
577 if ((ptr = strchr(line,'#'))) *ptr = '\0';
578 nwords = utils::count_words(line);
579 }
580
581 if (nwords != params_per_line)
582 error->all(FLERR,"Insufficient spline parameters in potential file");
583
584 // words = ptrs to all words in line
585
586 nwords = 0;
587 words[nwords++] = strtok(line," \t\n\r\f");
588 while ((words[nwords++] = strtok(nullptr," \t\n\r\f"))) continue;
589
590 // ielement,jelement,kelement = 1st args
591 // if all 3 args are in element list, then parse this line
592 // else skip to next line
593
594 for (ielement = 0; ielement < nelements; ielement++)
595 if (strcmp(words[0],elements[ielement]) == 0) break;
596 if (ielement == nelements) continue;
597 for (jelement = 0; jelement < nelements; jelement++)
598 if (strcmp(words[1],elements[jelement]) == 0) break;
599 if (jelement == nelements) continue;
600 for (kelement = 0; kelement < nelements; kelement++)
601 if (strcmp(words[2],elements[kelement]) == 0) break;
602 if (kelement == nelements) continue;
603
604 // load up parameter settings and error check their values
605
606 if (nparams == maxparam) {
607 maxparam += DELTA;
608 params = (Param *) memory->srealloc(params,maxparam*sizeof(Param),
609 "pair:params");
610
611 // make certain all addional allocated storage is initialized
612 // to avoid false positives when checking with valgrind
613
614 memset(params + nparams, 0, DELTA*sizeof(Param));
615 }
616
617 params[nparams].ielement = ielement;
618 params[nparams].jelement = jelement;
619 params[nparams].kelement = kelement;
620 params[nparams].powerm = atof(words[3]);
621 params[nparams].gamma = atof(words[4]);
622 params[nparams].lam3 = atof(words[5]);
623 params[nparams].c = atof(words[6]);
624 params[nparams].d = atof(words[7]);
625 params[nparams].h = atof(words[8]);
626 params[nparams].powern = atof(words[9]);
627 params[nparams].beta = atof(words[10]);
628 params[nparams].lam2 = atof(words[11]);
629 params[nparams].bigb = atof(words[12]);
630 params[nparams].bigr = atof(words[13]);
631 params[nparams].bigd = atof(words[14]);
632 params[nparams].lam1 = atof(words[15]);
633 params[nparams].biga = atof(words[16]);
634
635 // currently only allow m exponent of 1 or 3
636
637 params[nparams].powermint = int(params[nparams].powerm);
638
639 if (params[nparams].c < 0.0 || params[nparams].d < 0.0 ||
640 params[nparams].powern < 0.0 || params[nparams].beta < 0.0 ||
641 params[nparams].lam2 < 0.0 || params[nparams].bigb < 0.0 ||
642 params[nparams].bigr < 0.0 ||params[nparams].bigd < 0.0 ||
643 params[nparams].bigd > params[nparams].bigr ||
644 params[nparams].lam1 < 0.0 || params[nparams].biga < 0.0 ||
645 params[nparams].powerm - params[nparams].powermint != 0.0 ||
646 (params[nparams].powermint != 3 && params[nparams].powermint != 1) ||
647 params[nparams].gamma < 0.0)
648 error->all(FLERR,"Illegal ExTeP parameter");
649
650 nparams++;
651 if (nparams >= pow(nelements,3)) break;
652 }
653
654 // deallocate words array
655 delete [] words;
656
657 /* F_IJ (3) */
658 // read the spline coefficients
659 params_per_line = 8;
660 // reallocate with new size
661 words = new char*[params_per_line+1];
662
663 // initialize F_corr_data to all zeros
664 for (int iel=0;iel<nelements;iel++)
665 for (int jel=0;jel<nelements;jel++)
666 for (int in=0;in<4;in++)
667 for (int jn=0;jn<4;jn++)
668 for (int ivar=0;ivar<3;ivar++)
669 F_corr_data[iel][jel][in][jn][ivar]=0;
670
671 // loop until EOF
672 while (1) {
673 if (comm->me == 0) {
674 ptr = fgets(line,MAXLINE,fp);
675 //fputs(line,stdout);
676 if (ptr == nullptr) {
677 eof = 1;
678 fclose(fp);
679 } else n = strlen(line) + 1;
680 }
681 MPI_Bcast(&eof,1,MPI_INT,0,world);
682 if (eof) break;
683 MPI_Bcast(&n,1,MPI_INT,0,world);
684 MPI_Bcast(line,n,MPI_CHAR,0,world);
685
686 // strip comment, skip line if blank
687
688 if ((ptr = strchr(line,'#'))) *ptr = '\0';
689 nwords = utils::count_words(line);
690 if (nwords == 0) continue;
691
692 // words = ptrs to all words in line
693
694 nwords = 0;
695 words[nwords++] = strtok(line," \t\n\r\f");
696 while ((nwords < params_per_line)
697 && (words[nwords++] = strtok(nullptr," \t\n\r\f"))) continue;
698
699 // skip line if it is a leftover from the previous section,
700 // which can be identified by having 3 elements (instead of 2)
701 // as first words.
702
703 if (isupper(words[0][0]) && isupper(words[1][0]) && isupper(words[2][0]))
704 continue;
705
706 // need to have two elements followed by a number in each line
707 if (!(isupper(words[0][0]) && isupper(words[1][0])
708 && !isupper(words[2][0])))
709 error->all(FLERR,"Incorrect format in ExTeP potential file");
710
711 // ielement,jelement = 1st args
712 // if all 3 args are in element list, then parse this line
713 // else skip to next line
714 // these lines set ielement and jelement to the
715 // integers matching the strings from the input
716
717 for (ielement = 0; ielement < nelements; ielement++)
718 if (strcmp(words[0],elements[ielement]) == 0) break;
719 if (ielement == nelements) continue;
720 for (jelement = 0; jelement < nelements; jelement++)
721 if (strcmp(words[1],elements[jelement]) == 0) break;
722 if (jelement == nelements) continue;
723
724 int Ni = atoi(words[2]);
725 int Nj = atoi(words[3]);
726 double spline_val = atof(words[4]);
727 double spline_derx = atof(words[5]);
728 double spline_dery = atof(words[6]);
729
730 // Set value for all pairs of ielement,jelement (any kelement)
731 for (int iparam = 0; iparam < nparams; iparam++) {
732 if ( ielement == params[iparam].ielement
733 && jelement == params[iparam].jelement) {
734 F_corr_data[ielement][jelement][Ni][Nj][0] = spline_val;
735 F_corr_data[ielement][jelement][Ni][Nj][1] = spline_derx;
736 F_corr_data[ielement][jelement][Ni][Nj][2] = spline_dery;
737
738 F_corr_data[jelement][ielement][Nj][Ni][0] = spline_val;
739 F_corr_data[jelement][ielement][Nj][Ni][1] = spline_dery;
740 F_corr_data[jelement][ielement][Nj][Ni][2] = spline_derx;
741 }
742 }
743 }
744
745 delete [] words;
746 /* END F_IJ (3) */
747
748 }
749
750 /* ---------------------------------------------------------------------- */
751
setup()752 void PairExTeP::setup()
753 {
754 int i,j,k,m,n;
755
756 // set elem3param for all element triplet combinations
757 // must be a single exact match to lines read from file
758 // do not allow for ACB in place of ABC
759
760 memory->destroy(elem3param);
761 memory->create(elem3param,nelements,nelements,nelements,"pair:elem3param");
762
763 for (i = 0; i < nelements; i++)
764 for (j = 0; j < nelements; j++)
765 for (k = 0; k < nelements; k++) {
766 n = -1;
767 for (m = 0; m < nparams; m++) {
768 if (i == params[m].ielement && j == params[m].jelement &&
769 k == params[m].kelement) {
770 if (n >= 0) error->all(FLERR,"Potential file has duplicate entry");
771 n = m;
772 }
773 }
774 if (n < 0) error->all(FLERR,"Potential file is missing an entry");
775 elem3param[i][j][k] = n;
776 }
777
778 // compute parameter values derived from inputs
779
780 for (m = 0; m < nparams; m++) {
781 params[m].cut = params[m].bigr + params[m].bigd;
782 params[m].cutsq = params[m].cut*params[m].cut;
783
784 params[m].c1 = pow(2.0*params[m].powern*1.0e-16,-1.0/params[m].powern);
785 params[m].c2 = pow(2.0*params[m].powern*1.0e-8,-1.0/params[m].powern);
786 params[m].c3 = 1.0/params[m].c2;
787 params[m].c4 = 1.0/params[m].c1;
788 }
789
790 // set cutmax to max of all params
791
792 cutmax = 0.0;
793 for (m = 0; m < nparams; m++)
794 if (params[m].cut > cutmax) cutmax = params[m].cut;
795 }
796
797 /* ---------------------------------------------------------------------- */
798
repulsive(Param * param,double rsq,double & fforce,int eflag,double & eng)799 void PairExTeP::repulsive(Param *param, double rsq, double &fforce,
800 int eflag, double &eng)
801 {
802 double r,tmp_fc,tmp_fc_d,tmp_exp;
803
804 r = sqrt(rsq);
805 tmp_fc = ters_fc(r,param);
806 tmp_fc_d = ters_fc_d(r,param);
807 tmp_exp = exp(-param->lam1 * r);
808 fforce = -param->biga * tmp_exp * (tmp_fc_d - tmp_fc*param->lam1) / r;
809 if (eflag) eng = tmp_fc * param->biga * tmp_exp;
810 }
811
812 /* ---------------------------------------------------------------------- */
813
zeta(Param * param,double rij,double rik,double * delrij,double * delrik)814 double PairExTeP::zeta(Param *param, double rij, double rik,
815 double *delrij, double *delrik)
816 {
817 double costheta,arg,ex_delr;
818
819 costheta = (delrij[0]*delrik[0] + delrij[1]*delrik[1] +
820 delrij[2]*delrik[2]) / (rij*rik);
821
822 if (param->powermint == 3) arg = pow(param->lam3 * (rij-rik),3.0);
823 else arg = param->lam3 * (rij-rik);
824
825 if (arg > 69.0776) ex_delr = 1.e30;
826 else if (arg < -69.0776) ex_delr = 0.0;
827 else ex_delr = exp(arg);
828
829 return ters_fc(rik,param) * ters_gijk(costheta,param) * ex_delr;
830 }
831
832 /* ---------------------------------------------------------------------- */
833
force_zeta(Param * param,double r,double zeta_ij,double & fforce,double & prefactor,int eflag,double & eng)834 void PairExTeP::force_zeta(Param *param, double r, double zeta_ij,
835 double &fforce, double &prefactor,
836 int eflag, double &eng)
837 {
838 double fa,fa_d,bij;
839
840 fa = ters_fa(r,param);
841 fa_d = ters_fa_d(r,param);
842 bij = ters_bij(zeta_ij,param);
843 fforce = 0.5*bij*fa_d / r;
844 prefactor = -0.5*fa * ( ters_bij_d(zeta_ij,param) );
845 if (eflag) eng = 0.5*bij*fa;
846 }
847
848 /* ----------------------------------------------------------------------
849 attractive term
850 use param_ij cutoff for rij test
851 use param_ijk cutoff for rik test
852 ------------------------------------------------------------------------- */
853
attractive(Param * param,double prefactor,double rsqij,double rsqik,double * delrij,double * delrik,double * fi,double * fj,double * fk)854 void PairExTeP::attractive(Param *param, double prefactor,
855 double rsqij, double rsqik,
856 double *delrij, double *delrik,
857 double *fi, double *fj, double *fk)
858 {
859 double rij_hat[3],rik_hat[3];
860 double rij,rijinv,rik,rikinv;
861
862 rij = sqrt(rsqij);
863 rijinv = 1.0/rij;
864 scale3(rijinv,delrij,rij_hat);
865
866 rik = sqrt(rsqik);
867 rikinv = 1.0/rik;
868 scale3(rikinv,delrik,rik_hat);
869
870 ters_zetaterm_d(prefactor,rij_hat,rij,rik_hat,rik,fi,fj,fk,param);
871 }
872
873 /* ---------------------------------------------------------------------- */
874
ters_fc(double r,Param * param)875 double PairExTeP::ters_fc(double r, Param *param)
876 {
877 double ters_R = param->bigr;
878 double ters_D = param->bigd;
879
880 if (r < ters_R-ters_D) return 1.0;
881 if (r > ters_R+ters_D) return 0.0;
882 return 0.5*(1.0 - sin(MY_PI2*(r - ters_R)/ters_D));
883 }
884
885 /* ---------------------------------------------------------------------- */
886
ters_fc_d(double r,Param * param)887 double PairExTeP::ters_fc_d(double r, Param *param)
888 {
889 double ters_R = param->bigr;
890 double ters_D = param->bigd;
891
892 if (r < ters_R-ters_D) return 0.0;
893 if (r > ters_R+ters_D) return 0.0;
894 return -(MY_PI4/ters_D) * cos(MY_PI2*(r - ters_R)/ters_D);
895 }
896
897 /* ---------------------------------------------------------------------- */
898
ters_fa(double r,Param * param)899 double PairExTeP::ters_fa(double r, Param *param)
900 {
901 if (r > param->bigr + param->bigd) return 0.0;
902 return -param->bigb * exp(-param->lam2 * r) * ters_fc(r,param);
903 }
904
905 /* ---------------------------------------------------------------------- */
906
ters_fa_d(double r,Param * param)907 double PairExTeP::ters_fa_d(double r, Param *param)
908 {
909 if (r > param->bigr + param->bigd) return 0.0;
910 return param->bigb * exp(-param->lam2 * r) *
911 (param->lam2 * ters_fc(r,param) - ters_fc_d(r,param));
912 }
913
914 /* ---------------------------------------------------------------------- */
915
ters_bij(double zeta,Param * param)916 double PairExTeP::ters_bij(double zeta, Param *param)
917 {
918 double tmp = param->beta * zeta;
919 if (tmp > param->c1) return 1.0/sqrt(tmp);
920 if (tmp > param->c2)
921 return (1.0 - pow(tmp,-param->powern) / (2.0*param->powern))/sqrt(tmp);
922 if (tmp < param->c4) return 1.0;
923 if (tmp < param->c3)
924 return 1.0 - pow(tmp,param->powern)/(2.0*param->powern);
925 return pow(1.0 + pow(tmp,param->powern), -1.0/(2.0*param->powern));
926 }
927
928 /* ---------------------------------------------------------------------- */
929
ters_bij_d(double zeta,Param * param)930 double PairExTeP::ters_bij_d(double zeta, Param *param)
931 {
932 double tmp = param->beta * zeta;
933 if (tmp > param->c1) return param->beta * -0.5*pow(tmp,-1.5);
934 if (tmp > param->c2)
935 return param->beta * (-0.5*pow(tmp,-1.5) *
936 (1.0 - 0.5*(1.0 + 1.0/(2.0*param->powern)) *
937 pow(tmp,-param->powern)));
938 if (tmp < param->c4) return 0.0;
939 if (tmp < param->c3)
940 return -0.5*param->beta * pow(tmp,param->powern-1.0);
941
942 double tmp_n = pow(tmp,param->powern);
943 return -0.5 * pow(1.0+tmp_n, -1.0-(1.0/(2.0*param->powern)))*tmp_n / zeta;
944 }
945
946 /* ---------------------------------------------------------------------- */
947
ters_zetaterm_d(double prefactor,double * rij_hat,double rij,double * rik_hat,double rik,double * dri,double * drj,double * drk,Param * param)948 void PairExTeP::ters_zetaterm_d(double prefactor,
949 double *rij_hat, double rij,
950 double *rik_hat, double rik,
951 double *dri, double *drj, double *drk,
952 Param *param)
953 {
954 double gijk,gijk_d,ex_delr,ex_delr_d,fc,dfc,cos_theta,tmp;
955 double dcosdri[3],dcosdrj[3],dcosdrk[3];
956
957 fc = ters_fc(rik,param);
958 dfc = ters_fc_d(rik,param);
959 if (param->powermint == 3) tmp = pow(param->lam3 * (rij-rik),3.0);
960 else tmp = param->lam3 * (rij-rik);
961
962 if (tmp > 69.0776) ex_delr = 1.e30;
963 else if (tmp < -69.0776) ex_delr = 0.0;
964 else ex_delr = exp(tmp);
965
966 if (param->powermint == 3)
967 ex_delr_d = 3.0*pow(param->lam3,3.0) * pow(rij-rik,2.0)*ex_delr;
968 else ex_delr_d = param->lam3 * ex_delr;
969
970 cos_theta = dot3(rij_hat,rik_hat);
971 gijk = ters_gijk(cos_theta,param);
972 gijk_d = ters_gijk_d(cos_theta,param);
973 costheta_d(rij_hat,rij,rik_hat,rik,dcosdri,dcosdrj,dcosdrk);
974
975 // compute the derivative wrt Ri
976 // dri = -dfc*gijk*ex_delr*rik_hat;
977 // dri += fc*gijk_d*ex_delr*dcosdri;
978 // dri += fc*gijk*ex_delr_d*(rik_hat - rij_hat);
979
980 scale3(-dfc*gijk*ex_delr,rik_hat,dri);
981 scaleadd3(fc*gijk_d*ex_delr,dcosdri,dri,dri);
982 scaleadd3(fc*gijk*ex_delr_d,rik_hat,dri,dri);
983 scaleadd3(-fc*gijk*ex_delr_d,rij_hat,dri,dri);
984 scale3(prefactor,dri);
985
986 // compute the derivative wrt Rj
987 // drj = fc*gijk_d*ex_delr*dcosdrj;
988 // drj += fc*gijk*ex_delr_d*rij_hat;
989
990 scale3(fc*gijk_d*ex_delr,dcosdrj,drj);
991 scaleadd3(fc*gijk*ex_delr_d,rij_hat,drj,drj);
992 scale3(prefactor,drj);
993
994 // compute the derivative wrt Rk
995 // drk = dfc*gijk*ex_delr*rik_hat;
996 // drk += fc*gijk_d*ex_delr*dcosdrk;
997 // drk += -fc*gijk*ex_delr_d*rik_hat;
998
999 scale3(dfc*gijk*ex_delr,rik_hat,drk);
1000 scaleadd3(fc*gijk_d*ex_delr,dcosdrk,drk,drk);
1001 scaleadd3(-fc*gijk*ex_delr_d,rik_hat,drk,drk);
1002 scale3(prefactor,drk);
1003 }
1004
1005 /* ---------------------------------------------------------------------- */
1006
costheta_d(double * rij_hat,double rij,double * rik_hat,double rik,double * dri,double * drj,double * drk)1007 void PairExTeP::costheta_d(double *rij_hat, double rij,
1008 double *rik_hat, double rik,
1009 double *dri, double *drj, double *drk)
1010 {
1011 // first element is devative wrt Ri, second wrt Rj, third wrt Rk
1012
1013 double cos_theta = dot3(rij_hat,rik_hat);
1014
1015 scaleadd3(-cos_theta,rij_hat,rik_hat,drj);
1016 scale3(1.0/rij,drj);
1017 scaleadd3(-cos_theta,rik_hat,rij_hat,drk);
1018 scale3(1.0/rik,drk);
1019 add3(drj,drk,dri);
1020 scale3(-1.0,dri);
1021 }
1022
1023
1024 /* ---------------------------------------------------------------------- */
1025
1026 /* F_IJ (4) */
1027 // initialize spline for F_corr (based on PairLCBOP::F_conj)
1028
spline_init()1029 void PairExTeP::spline_init() {
1030 for ( int iel=0; iel<nelements; iel++) {
1031 for ( int jel=0; jel<nelements; jel++) {
1032 for (int N_ij=0; N_ij<4; N_ij++) {
1033 for (int N_ji=0; N_ji<4; N_ji++) {
1034 TF_corr_param &f = F_corr_param[iel][jel][N_ij][N_ji];
1035
1036 // corner points for each spline function
1037 f.f_00 = F_corr_data[iel][jel][N_ij ][N_ji ][0];
1038 f.f_01 = F_corr_data[iel][jel][N_ij ][N_ji+1][0];
1039 f.f_10 = F_corr_data[iel][jel][N_ij+1][N_ji ][0];
1040 f.f_11 = F_corr_data[iel][jel][N_ij+1][N_ji+1][0];
1041
1042 // compute f tilde according to [Los & Fasolino, PRB 68, 024107 2003]
1043 f.f_x_00 = F_corr_data[iel][jel][N_ij ][N_ji ][1] - f.f_10 + f.f_00;
1044 f.f_x_01 = F_corr_data[iel][jel][N_ij ][N_ji+1][1] - f.f_11 + f.f_01;
1045 f.f_x_10 = -(F_corr_data[iel][jel][N_ij+1][N_ji ][1] - f.f_10 + f.f_00);
1046 f.f_x_11 = -(F_corr_data[iel][jel][N_ij+1][N_ji+1][1] - f.f_11 + f.f_01);
1047
1048 f.f_y_00 = F_corr_data[iel][jel][N_ij ][N_ji ][2] - f.f_01 + f.f_00;
1049 f.f_y_01 = -(F_corr_data[iel][jel][N_ij ][N_ji+1][2] - f.f_01 + f.f_00);
1050 f.f_y_10 = F_corr_data[iel][jel][N_ij+1][N_ji ][2] - f.f_11 + f.f_10;
1051 f.f_y_11 = -(F_corr_data[iel][jel][N_ij+1][N_ji+1][2] - f.f_11 + f.f_10);
1052 }
1053 }
1054 }
1055 }
1056 }
1057
envelop_function(double x,double y,double * func_der)1058 double PairExTeP::envelop_function(double x, double y, double *func_der) {
1059 double fx,fy,fxy,dfx,dfxydx;
1060 double del, delsq;
1061
1062 fxy = 1.0;
1063 dfxydx = 0.0;
1064
1065 if (x <= 3.0) {
1066 fx = 1.0;
1067 dfx = 0.0;
1068 if (x < 1.0 && y < 1.0) {
1069 double gx=(1.0-x);
1070 double gy=(1.0-y);
1071 double gxsq=gx*gx;
1072 double gysq=gy*gy;
1073 fxy = 1.0 - gxsq*gysq;
1074 dfxydx = 2.0*gx*gysq;
1075 }
1076 } else if (x < 4.0) {
1077 del = 4.0-x;
1078 delsq = del*del;
1079 fx = (3.0-2.0*del)*delsq;
1080 dfx = - 6.0*del*(1.0-del);
1081 } else {
1082 fx = 0.0;
1083 dfx = 0.0;
1084 }
1085 if (y <= 3.0) {
1086 fy = 1.0;
1087 } else if (y < 4.0) {
1088 del = 4.0-y;
1089 delsq = del*del;
1090 fy = (3.0-2.0*del)*delsq;
1091 } else {
1092 fy = 0.0;
1093 }
1094
1095 double func_val = fxy*fx*fy;
1096 *func_der = dfxydx*fx*fy+fxy*dfx*fy;
1097
1098 return func_val;
1099 }
1100
F_corr(int iel,int jel,double Ndij,double Ndji,double * dFN_x,double * dFN_y)1101 double PairExTeP::F_corr(int iel, int jel, double Ndij, double Ndji, double *dFN_x, double *dFN_y) {
1102
1103 // compute F_XY
1104
1105 int Ndij_int = static_cast<int>( floor( Ndij ) );
1106 int Ndji_int = static_cast<int>( floor( Ndji ) );
1107 double x = Ndij - Ndij_int;
1108 double y = Ndji - Ndji_int;
1109 TF_corr_param &f = F_corr_param[iel][jel][Ndij_int][Ndji_int];
1110 double F = 0;
1111 double dF_dx = 0, dF_dy = 0;
1112 double l, r;
1113 if (Ndij_int < 4 && Ndji_int < 4) {
1114 l = (1-y)* (1-x);
1115 r = ( f.f_00 + x*x* f.f_x_10 + y*y* f.f_y_01 );
1116 F += l*r;
1117 dF_dx += -(1-y)*r +l*2*x* f.f_x_10;
1118 dF_dy += -(1-x)*r +l*2*y* f.f_y_01;
1119 l = (1-y)*x;
1120 r = ( f.f_10 + (1-x)*(1-x)*f.f_x_00 + y* y* f.f_y_11 );
1121 F += l*r;
1122 dF_dx += (1-y)*r -l*2*(1-x)*f.f_x_00;
1123 dF_dy += -x*r +l*2*y* f.f_y_11;
1124 l = y* (1-x);
1125 r = ( f.f_01 + x*x* f.f_x_11 + (1-y)*(1-y)*f.f_y_00 );
1126 F += l*r;
1127 dF_dx += -y*r +l*2*x* f.f_x_11;
1128 dF_dy += (1-x)*r -l*2*(1-y)*f.f_y_00;
1129 l = y* x;
1130 r = ( f.f_11 + (1-x)*(1-x)*f.f_x_01 + (1-y)*(1-y)*f.f_y_10 );
1131 F += l*r;
1132 dF_dx += y*r -l*2*(1-x)*f.f_x_01;
1133 dF_dy += x*r -l*2*(1-y)*f.f_y_10;
1134 }
1135 double result = F;
1136 *dFN_x = dF_dx;
1137 *dFN_y = dF_dy;
1138
1139 return result;
1140 }
1141 /* F_IJ (4) */
1142