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: Paul Crozier (SNL)
17 The lj-fsw/coul-fsh (force-switched and force-shifted) sections
18 were provided by Robert Meissner
19 and Lucio Colombi Ciacchi of Bremen University, Bremen, Germany,
20 with additional assistance from Robert A. Latour, Clemson University
21 ------------------------------------------------------------------------- */
22
23 #include "pair_lj_charmmfsw_coul_charmmfsh.h"
24
25 #include <cmath>
26 #include <cstring>
27 #include "atom.h"
28 #include "update.h"
29 #include "comm.h"
30 #include "force.h"
31 #include "neighbor.h"
32 #include "neigh_list.h"
33 #include "memory.h"
34 #include "error.h"
35
36
37 using namespace LAMMPS_NS;
38
39 /* ---------------------------------------------------------------------- */
40
PairLJCharmmfswCoulCharmmfsh(LAMMPS * lmp)41 PairLJCharmmfswCoulCharmmfsh::PairLJCharmmfswCoulCharmmfsh(LAMMPS *lmp) :
42 Pair(lmp)
43 {
44 implicit = 0;
45 mix_flag = ARITHMETIC;
46 writedata = 1;
47
48 // short-range/long-range flag accessed by DihedralCharmmfsw
49
50 dihedflag = 0;
51
52 // switch qqr2e from LAMMPS value to CHARMM value
53
54 if (strcmp(update->unit_style,"real") == 0) {
55 if ((comm->me == 0) && (force->qqr2e != force->qqr2e_charmm_real))
56 error->message(FLERR,"Switching to CHARMM coulomb energy"
57 " conversion constant");
58 force->qqr2e = force->qqr2e_charmm_real;
59 }
60 }
61
62 /* ---------------------------------------------------------------------- */
63
~PairLJCharmmfswCoulCharmmfsh()64 PairLJCharmmfswCoulCharmmfsh::~PairLJCharmmfswCoulCharmmfsh()
65 {
66 // switch qqr2e back from CHARMM value to LAMMPS value
67
68 if (update && strcmp(update->unit_style,"real") == 0) {
69 if ((comm->me == 0) && (force->qqr2e == force->qqr2e_charmm_real))
70 error->message(FLERR,"Restoring original LAMMPS coulomb energy"
71 " conversion constant");
72 force->qqr2e = force->qqr2e_lammps_real;
73 }
74
75 if (copymode) return;
76
77 if (allocated) {
78 memory->destroy(setflag);
79 memory->destroy(cutsq);
80
81 memory->destroy(epsilon);
82 memory->destroy(sigma);
83 memory->destroy(eps14);
84 memory->destroy(sigma14);
85 memory->destroy(lj1);
86 memory->destroy(lj2);
87 memory->destroy(lj3);
88 memory->destroy(lj4);
89 memory->destroy(lj14_1);
90 memory->destroy(lj14_2);
91 memory->destroy(lj14_3);
92 memory->destroy(lj14_4);
93 }
94 }
95
96 /* ---------------------------------------------------------------------- */
97
compute(int eflag,int vflag)98 void PairLJCharmmfswCoulCharmmfsh::compute(int eflag, int vflag)
99 {
100 int i,j,ii,jj,inum,jnum,itype,jtype;
101 double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,evdwl,evdwl12,evdwl6,ecoul,fpair;
102 double r,rinv,r3inv,rsq,r2inv,r6inv,forcecoul,forcelj,factor_coul,factor_lj;
103 double switch1;
104 int *ilist,*jlist,*numneigh,**firstneigh;
105
106 evdwl = ecoul = 0.0;
107 ev_init(eflag,vflag);
108
109 double **x = atom->x;
110 double **f = atom->f;
111 double *q = atom->q;
112 int *type = atom->type;
113 int nlocal = atom->nlocal;
114 double *special_coul = force->special_coul;
115 double *special_lj = force->special_lj;
116 int newton_pair = force->newton_pair;
117 double qqrd2e = force->qqrd2e;
118
119 inum = list->inum;
120 ilist = list->ilist;
121 numneigh = list->numneigh;
122 firstneigh = list->firstneigh;
123
124 // loop over neighbors of my atoms
125
126 for (ii = 0; ii < inum; ii++) {
127 i = ilist[ii];
128 qtmp = q[i];
129 xtmp = x[i][0];
130 ytmp = x[i][1];
131 ztmp = x[i][2];
132 itype = type[i];
133 jlist = firstneigh[i];
134 jnum = numneigh[i];
135
136 for (jj = 0; jj < jnum; jj++) {
137 j = jlist[jj];
138 factor_lj = special_lj[sbmask(j)];
139 factor_coul = special_coul[sbmask(j)];
140 j &= NEIGHMASK;
141
142 delx = xtmp - x[j][0];
143 dely = ytmp - x[j][1];
144 delz = ztmp - x[j][2];
145 rsq = delx*delx + dely*dely + delz*delz;
146
147 if (rsq < cut_bothsq) {
148 r2inv = 1.0/rsq;
149 r = sqrt(rsq);
150
151 if (rsq < cut_coulsq) {
152 forcecoul = qqrd2e * qtmp*q[j]*
153 (sqrt(r2inv) - r*cut_coulinv*cut_coulinv);
154 } else forcecoul = 0.0;
155
156 if (rsq < cut_ljsq) {
157 r6inv = r2inv*r2inv*r2inv;
158 jtype = type[j];
159 forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]);
160 if (rsq > cut_lj_innersq) {
161 switch1 = (cut_ljsq-rsq) * (cut_ljsq-rsq) *
162 (cut_ljsq + 2.0*rsq - 3.0*cut_lj_innersq) / denom_lj;
163 forcelj = forcelj*switch1;
164 }
165 } else forcelj = 0.0;
166
167 fpair = (factor_coul*forcecoul + factor_lj*forcelj) * r2inv;
168
169 f[i][0] += delx*fpair;
170 f[i][1] += dely*fpair;
171 f[i][2] += delz*fpair;
172 if (newton_pair || j < nlocal) {
173 f[j][0] -= delx*fpair;
174 f[j][1] -= dely*fpair;
175 f[j][2] -= delz*fpair;
176 }
177
178 if (eflag) {
179 if (rsq < cut_coulsq) {
180 ecoul = qqrd2e * qtmp*q[j]*
181 (sqrt(r2inv) + cut_coulinv*cut_coulinv*r - 2.0*cut_coulinv);
182 ecoul *= factor_coul;
183 } else ecoul = 0.0;
184 if (rsq < cut_ljsq) {
185 if (rsq > cut_lj_innersq) {
186 rinv = 1.0/r;
187 r3inv = rinv*rinv*rinv;
188 evdwl12 = lj3[itype][jtype]*cut_lj6*denom_lj12 *
189 (r6inv - cut_lj6inv)*(r6inv - cut_lj6inv);
190 evdwl6 = -lj4[itype][jtype]*cut_lj3*denom_lj6 *
191 (r3inv - cut_lj3inv)*(r3inv - cut_lj3inv);;
192 evdwl = evdwl12 + evdwl6;
193 } else {
194 evdwl12 = r6inv*lj3[itype][jtype]*r6inv -
195 lj3[itype][jtype]*cut_lj_inner6inv*cut_lj6inv;
196 evdwl6 = -lj4[itype][jtype]*r6inv +
197 lj4[itype][jtype]*cut_lj_inner3inv*cut_lj3inv;
198 evdwl = evdwl12 + evdwl6;
199 }
200 evdwl *= factor_lj;
201 } else evdwl = 0.0;
202 }
203
204 if (evflag) ev_tally(i,j,nlocal,newton_pair,
205 evdwl,ecoul,fpair,delx,dely,delz);
206 }
207 }
208 }
209
210 if (vflag_fdotr) virial_fdotr_compute();
211 }
212
213 /* ----------------------------------------------------------------------
214 allocate all arrays
215 ------------------------------------------------------------------------- */
216
allocate()217 void PairLJCharmmfswCoulCharmmfsh::allocate()
218 {
219 allocated = 1;
220 int n = atom->ntypes;
221
222 memory->create(setflag,n+1,n+1,"pair:setflag");
223 for (int i = 1; i <= n; i++)
224 for (int j = i; j <= n; j++)
225 setflag[i][j] = 0;
226
227 memory->create(cutsq,n+1,n+1,"pair:cutsq");
228
229 memory->create(epsilon,n+1,n+1,"pair:epsilon");
230 memory->create(sigma,n+1,n+1,"pair:sigma");
231 memory->create(eps14,n+1,n+1,"pair:eps14");
232 memory->create(sigma14,n+1,n+1,"pair:sigma14");
233 memory->create(lj1,n+1,n+1,"pair:lj1");
234 memory->create(lj2,n+1,n+1,"pair:lj2");
235 memory->create(lj3,n+1,n+1,"pair:lj3");
236 memory->create(lj4,n+1,n+1,"pair:lj4");
237 memory->create(lj14_1,n+1,n+1,"pair:lj14_1");
238 memory->create(lj14_2,n+1,n+1,"pair:lj14_2");
239 memory->create(lj14_3,n+1,n+1,"pair:lj14_3");
240 memory->create(lj14_4,n+1,n+1,"pair:lj14_4");
241 }
242
243 /* ----------------------------------------------------------------------
244 global settings
245 unlike other pair styles,
246 there are no individual pair settings that these override
247 ------------------------------------------------------------------------- */
248
settings(int narg,char ** arg)249 void PairLJCharmmfswCoulCharmmfsh::settings(int narg, char **arg)
250 {
251 if (narg != 2 && narg != 3)
252 error->all(FLERR,"Illegal pair_style command");
253
254 cut_lj_inner = utils::numeric(FLERR,arg[0],false,lmp);
255 cut_lj = utils::numeric(FLERR,arg[1],false,lmp);
256 if (narg == 2) {
257 cut_coul = cut_lj;
258 } else {
259 cut_coul = utils::numeric(FLERR,arg[2],false,lmp);
260 }
261 }
262
263 /* ----------------------------------------------------------------------
264 set coeffs for one or more type pairs
265 ------------------------------------------------------------------------- */
266
coeff(int narg,char ** arg)267 void PairLJCharmmfswCoulCharmmfsh::coeff(int narg, char **arg)
268 {
269 if (narg != 4 && narg != 6)
270 error->all(FLERR,"Incorrect args for pair coefficients");
271 if (!allocated) allocate();
272
273 int ilo,ihi,jlo,jhi;
274 utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error);
275 utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error);
276
277 double epsilon_one = utils::numeric(FLERR,arg[2],false,lmp);
278 double sigma_one = utils::numeric(FLERR,arg[3],false,lmp);
279 double eps14_one = epsilon_one;
280 double sigma14_one = sigma_one;
281 if (narg == 6) {
282 eps14_one = utils::numeric(FLERR,arg[4],false,lmp);
283 sigma14_one = utils::numeric(FLERR,arg[5],false,lmp);
284 }
285
286 int count = 0;
287 for (int i = ilo; i <= ihi; i++) {
288 for (int j = MAX(jlo,i); j <= jhi; j++) {
289 epsilon[i][j] = epsilon_one;
290 sigma[i][j] = sigma_one;
291 eps14[i][j] = eps14_one;
292 sigma14[i][j] = sigma14_one;
293 setflag[i][j] = 1;
294 count++;
295 }
296 }
297
298 if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients");
299 }
300
301 /* ----------------------------------------------------------------------
302 init specific to this pair style
303 ------------------------------------------------------------------------- */
304
init_style()305 void PairLJCharmmfswCoulCharmmfsh::init_style()
306 {
307 if (!atom->q_flag)
308 error->all(FLERR,"Pair style lj/charmmfsw/coul/charmmfsh "
309 "requires atom attribute q");
310
311 neighbor->request(this,instance_me);
312
313 // require cut_lj_inner < cut_lj
314
315 if (cut_lj_inner >= cut_lj)
316 error->all(FLERR,"Pair inner lj cutoff >= Pair outer lj cutoff");
317
318 cut_lj_innersq = cut_lj_inner * cut_lj_inner;
319 cut_ljsq = cut_lj * cut_lj;
320 cut_ljinv = 1.0/cut_lj;
321 cut_lj_innerinv = 1.0/cut_lj_inner;
322 cut_lj3 = cut_lj * cut_lj * cut_lj;
323 cut_lj3inv = cut_ljinv * cut_ljinv * cut_ljinv;
324 cut_lj_inner3inv = cut_lj_innerinv * cut_lj_innerinv * cut_lj_innerinv;
325 cut_lj_inner3 = cut_lj_inner * cut_lj_inner * cut_lj_inner;
326 cut_lj6 = cut_ljsq * cut_ljsq * cut_ljsq;
327 cut_lj6inv = cut_lj3inv * cut_lj3inv;
328 cut_lj_inner6inv = cut_lj_inner3inv * cut_lj_inner3inv;
329 cut_lj_inner6 = cut_lj_innersq * cut_lj_innersq * cut_lj_innersq;
330 cut_coulsq = cut_coul * cut_coul;
331 cut_coulinv = 1.0/cut_coul;
332 cut_bothsq = MAX(cut_ljsq,cut_coulsq);
333
334 denom_lj = (cut_ljsq-cut_lj_innersq) * (cut_ljsq-cut_lj_innersq) *
335 (cut_ljsq-cut_lj_innersq);
336 denom_lj12 = 1.0/(cut_lj6 - cut_lj_inner6);
337 denom_lj6 = 1.0/(cut_lj3 - cut_lj_inner3);
338 }
339
340 /* ----------------------------------------------------------------------
341 init for one type pair i,j and corresponding j,i
342 ------------------------------------------------------------------------- */
343
init_one(int i,int j)344 double PairLJCharmmfswCoulCharmmfsh::init_one(int i, int j)
345 {
346 if (setflag[i][j] == 0) {
347 epsilon[i][j] = mix_energy(epsilon[i][i],epsilon[j][j],
348 sigma[i][i],sigma[j][j]);
349 sigma[i][j] = mix_distance(sigma[i][i],sigma[j][j]);
350 eps14[i][j] = mix_energy(eps14[i][i],eps14[j][j],
351 sigma14[i][i],sigma14[j][j]);
352 sigma14[i][j] = mix_distance(sigma14[i][i],sigma14[j][j]);
353 }
354
355 double cut = MAX(cut_lj,cut_coul);
356
357 lj1[i][j] = 48.0 * epsilon[i][j] * pow(sigma[i][j],12.0);
358 lj2[i][j] = 24.0 * epsilon[i][j] * pow(sigma[i][j],6.0);
359 lj3[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],12.0);
360 lj4[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],6.0);
361 lj14_1[i][j] = 48.0 * eps14[i][j] * pow(sigma14[i][j],12.0);
362 lj14_2[i][j] = 24.0 * eps14[i][j] * pow(sigma14[i][j],6.0);
363 lj14_3[i][j] = 4.0 * eps14[i][j] * pow(sigma14[i][j],12.0);
364 lj14_4[i][j] = 4.0 * eps14[i][j] * pow(sigma14[i][j],6.0);
365
366 lj1[j][i] = lj1[i][j];
367 lj2[j][i] = lj2[i][j];
368 lj3[j][i] = lj3[i][j];
369 lj4[j][i] = lj4[i][j];
370 lj14_1[j][i] = lj14_1[i][j];
371 lj14_2[j][i] = lj14_2[i][j];
372 lj14_3[j][i] = lj14_3[i][j];
373 lj14_4[j][i] = lj14_4[i][j];
374
375 return cut;
376 }
377
378 /* ----------------------------------------------------------------------
379 proc 0 writes to data file
380 ------------------------------------------------------------------------- */
381
write_data(FILE * fp)382 void PairLJCharmmfswCoulCharmmfsh::write_data(FILE *fp)
383 {
384 for (int i = 1; i <= atom->ntypes; i++)
385 fprintf(fp,"%d %g %g %g %g\n",
386 i,epsilon[i][i],sigma[i][i],eps14[i][i],sigma14[i][i]);
387 }
388
389 /* ----------------------------------------------------------------------
390 proc 0 writes all pairs to data file
391 ------------------------------------------------------------------------- */
392
write_data_all(FILE * fp)393 void PairLJCharmmfswCoulCharmmfsh::write_data_all(FILE *fp)
394 {
395 for (int i = 1; i <= atom->ntypes; i++)
396 for (int j = i; j <= atom->ntypes; j++)
397 fprintf(fp,"%d %d %g %g %g %g\n",i,j,
398 epsilon[i][j],sigma[i][j],eps14[i][j],sigma14[i][j]);
399 }
400
401
402 /* ----------------------------------------------------------------------
403 proc 0 writes to restart file
404 ------------------------------------------------------------------------- */
405
write_restart(FILE * fp)406 void PairLJCharmmfswCoulCharmmfsh::write_restart(FILE *fp)
407 {
408 write_restart_settings(fp);
409
410 int i,j;
411 for (i = 1; i <= atom->ntypes; i++)
412 for (j = i; j <= atom->ntypes; j++) {
413 fwrite(&setflag[i][j],sizeof(int),1,fp);
414 if (setflag[i][j]) {
415 fwrite(&epsilon[i][j],sizeof(double),1,fp);
416 fwrite(&sigma[i][j],sizeof(double),1,fp);
417 fwrite(&eps14[i][j],sizeof(double),1,fp);
418 fwrite(&sigma14[i][j],sizeof(double),1,fp);
419 }
420 }
421 }
422
423 /* ----------------------------------------------------------------------
424 proc 0 reads from restart file, bcasts
425 ------------------------------------------------------------------------- */
426
read_restart(FILE * fp)427 void PairLJCharmmfswCoulCharmmfsh::read_restart(FILE *fp)
428 {
429 read_restart_settings(fp);
430
431 allocate();
432
433 int i,j;
434 int me = comm->me;
435 for (i = 1; i <= atom->ntypes; i++)
436 for (j = i; j <= atom->ntypes; j++) {
437 if (me == 0) utils::sfread(FLERR,&setflag[i][j],sizeof(int),1,fp,nullptr,error);
438 MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world);
439 if (setflag[i][j]) {
440 if (me == 0) {
441 utils::sfread(FLERR,&epsilon[i][j],sizeof(double),1,fp,nullptr,error);
442 utils::sfread(FLERR,&sigma[i][j],sizeof(double),1,fp,nullptr,error);
443 utils::sfread(FLERR,&eps14[i][j],sizeof(double),1,fp,nullptr,error);
444 utils::sfread(FLERR,&sigma14[i][j],sizeof(double),1,fp,nullptr,error);
445 }
446 MPI_Bcast(&epsilon[i][j],1,MPI_DOUBLE,0,world);
447 MPI_Bcast(&sigma[i][j],1,MPI_DOUBLE,0,world);
448 MPI_Bcast(&eps14[i][j],1,MPI_DOUBLE,0,world);
449 MPI_Bcast(&sigma14[i][j],1,MPI_DOUBLE,0,world);
450 }
451 }
452 }
453
454 /* ----------------------------------------------------------------------
455 proc 0 writes to restart file
456 ------------------------------------------------------------------------- */
457
write_restart_settings(FILE * fp)458 void PairLJCharmmfswCoulCharmmfsh::write_restart_settings(FILE *fp)
459 {
460 fwrite(&cut_lj_inner,sizeof(double),1,fp);
461 fwrite(&cut_lj,sizeof(double),1,fp);
462 fwrite(&cut_coul,sizeof(double),1,fp);
463 fwrite(&offset_flag,sizeof(int),1,fp);
464 fwrite(&mix_flag,sizeof(int),1,fp);
465 }
466
467 /* ----------------------------------------------------------------------
468 proc 0 reads from restart file, bcasts
469 ------------------------------------------------------------------------- */
470
read_restart_settings(FILE * fp)471 void PairLJCharmmfswCoulCharmmfsh::read_restart_settings(FILE *fp)
472 {
473 if (comm->me == 0) {
474 utils::sfread(FLERR,&cut_lj_inner,sizeof(double),1,fp,nullptr,error);
475 utils::sfread(FLERR,&cut_lj,sizeof(double),1,fp,nullptr,error);
476 utils::sfread(FLERR,&cut_coul,sizeof(double),1,fp,nullptr,error);
477 utils::sfread(FLERR,&offset_flag,sizeof(int),1,fp,nullptr,error);
478 utils::sfread(FLERR,&mix_flag,sizeof(int),1,fp,nullptr,error);
479 }
480 MPI_Bcast(&cut_lj_inner,1,MPI_DOUBLE,0,world);
481 MPI_Bcast(&cut_lj,1,MPI_DOUBLE,0,world);
482 MPI_Bcast(&cut_coul,1,MPI_DOUBLE,0,world);
483 MPI_Bcast(&offset_flag,1,MPI_INT,0,world);
484 MPI_Bcast(&mix_flag,1,MPI_INT,0,world);
485 }
486
487 /* ---------------------------------------------------------------------- */
488
489 double PairLJCharmmfswCoulCharmmfsh::
single(int i,int j,int itype,int jtype,double rsq,double factor_coul,double factor_lj,double & fforce)490 single(int i, int j, int itype, int jtype,
491 double rsq, double factor_coul, double factor_lj, double &fforce)
492 {
493 double r,rinv,r2inv,r3inv,r6inv,forcecoul,forcelj;
494 double phicoul,philj,philj12,philj6;
495 double switch1;
496
497 r2inv = 1.0/rsq;
498 r = sqrt(rsq);
499 rinv = 1.0/r;
500 if (rsq < cut_coulsq) {
501 forcecoul = force->qqrd2e * atom->q[i]*atom->q[j] *
502 (sqrt(r2inv) - r*cut_coulinv*cut_coulinv);
503 } else forcecoul = 0.0;
504
505 if (rsq < cut_ljsq) {
506 r6inv = r2inv*r2inv*r2inv;
507 r3inv = rinv*rinv*rinv;
508 forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]);
509 if (rsq > cut_lj_innersq) {
510 switch1 = (cut_ljsq-rsq) * (cut_ljsq-rsq) *
511 (cut_ljsq + 2.0*rsq - 3.0*cut_lj_innersq) / denom_lj;
512 forcelj = forcelj*switch1;
513 }
514 } else forcelj = 0.0;
515
516 fforce = (factor_coul*forcecoul + factor_lj*forcelj) * r2inv;
517
518 double eng = 0.0;
519 if (rsq < cut_coulsq) {
520 phicoul = force->qqrd2e * atom->q[i]*atom->q[j] *
521 (sqrt(r2inv) + cut_coulinv*cut_coulinv*r - 2.0*cut_coulinv);
522 eng += factor_coul*phicoul;
523 }
524 if (rsq < cut_ljsq) {
525 if (rsq > cut_lj_innersq) {
526 philj12 = lj3[itype][jtype]*cut_lj6*denom_lj12 *
527 (r6inv - cut_lj6inv)*(r6inv - cut_lj6inv);
528 philj6 = -lj4[itype][jtype]*cut_lj3*denom_lj6 *
529 (r3inv - cut_lj3inv)*(r3inv - cut_lj3inv);;
530 philj = philj12 + philj6;
531 } else {
532 philj12 = r6inv*lj3[itype][jtype]*r6inv -
533 lj3[itype][jtype]*cut_lj_inner6inv*cut_lj6inv;
534 philj6 = -lj4[itype][jtype]*r6inv +
535 lj4[itype][jtype]*cut_lj_inner3inv*cut_lj3inv;
536 philj = philj12 + philj6;
537 }
538 eng += factor_lj*philj;
539 }
540
541 return eng;
542 }
543
544 /* ---------------------------------------------------------------------- */
545
extract(const char * str,int & dim)546 void *PairLJCharmmfswCoulCharmmfsh::extract(const char *str, int &dim)
547 {
548 dim = 2;
549 if (strcmp(str,"lj14_1") == 0) return (void *) lj14_1;
550 if (strcmp(str,"lj14_2") == 0) return (void *) lj14_2;
551 if (strcmp(str,"lj14_3") == 0) return (void *) lj14_3;
552 if (strcmp(str,"lj14_4") == 0) return (void *) lj14_4;
553
554 dim = 0;
555 if (strcmp(str,"implicit") == 0) return (void *) &implicit;
556
557 // info extracted by dihedral_charmmfsw
558
559 if (strcmp(str,"cut_coul") == 0) return (void *) &cut_coul;
560 if (strcmp(str,"cut_lj_inner") == 0) return (void *) &cut_lj_inner;
561 if (strcmp(str,"cut_lj") == 0) return (void *) &cut_lj;
562 if (strcmp(str,"dihedflag") == 0) return (void *) &dihedflag;
563
564 return nullptr;
565 }
566