lammps-sys 0.6.0

Generates bindings to LAMMPS' C interface (with optional builds from source)
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
/* ----------------------------------------------------------------------
   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.
------------------------------------------------------------------------- */

/* ----------------------------------------------------------------------
   Contributing author: Chuanfu Luo (luochuanfu@gmail.com)
------------------------------------------------------------------------- */

#include "bond_table.h"
#include <mpi.h>
#include <cmath>
#include <cstdlib>
#include <cstring>
#include "atom.h"
#include "neighbor.h"
#include "comm.h"
#include "force.h"
#include "memory.h"
#include "error.h"
#include "utils.h"

using namespace LAMMPS_NS;

enum{NONE,LINEAR,SPLINE};

#define MAXLINE 1024
#define BIGNUM 1.0e300

/* ---------------------------------------------------------------------- */

BondTable::BondTable(LAMMPS *lmp) : Bond(lmp)
{
  writedata = 0;
  ntables = 0;
  tables = NULL;
}

/* ---------------------------------------------------------------------- */

BondTable::~BondTable()
{
  for (int m = 0; m < ntables; m++) free_table(&tables[m]);
  memory->sfree(tables);

  if (allocated) {
    memory->destroy(setflag);
    memory->destroy(r0);
    memory->destroy(tabindex);
  }
}

/* ---------------------------------------------------------------------- */

void BondTable::compute(int eflag, int vflag)
{
  int i1,i2,n,type;
  double delx,dely,delz,ebond,fbond;
  double rsq,r;
  double u,mdu;

  ebond = 0.0;
  ev_init(eflag,vflag);

  double **x = atom->x;
  double **f = atom->f;
  int **bondlist = neighbor->bondlist;
  int nbondlist = neighbor->nbondlist;
  int nlocal = atom->nlocal;
  int newton_bond = force->newton_bond;

  for (n = 0; n < nbondlist; n++) {
    i1 = bondlist[n][0];
    i2 = bondlist[n][1];
    type = bondlist[n][2];

    delx = x[i1][0] - x[i2][0];
    dely = x[i1][1] - x[i2][1];
    delz = x[i1][2] - x[i2][2];

    rsq = delx*delx + dely*dely + delz*delz;
    r = sqrt(rsq);

    // force & energy

    uf_lookup(type,r,u,mdu);
    fbond = mdu/r;
    ebond = u;

    // apply force to each of 2 atoms

    if (newton_bond || i1 < nlocal) {
      f[i1][0] += delx*fbond;
      f[i1][1] += dely*fbond;
      f[i1][2] += delz*fbond;
    }

    if (newton_bond || i2 < nlocal) {
      f[i2][0] -= delx*fbond;
      f[i2][1] -= dely*fbond;
      f[i2][2] -= delz*fbond;
    }

    if (evflag) ev_tally(i1,i2,nlocal,newton_bond,ebond,fbond,delx,dely,delz);
  }
}

/* ---------------------------------------------------------------------- */

void BondTable::allocate()
{
  allocated = 1;
  int n = atom->nbondtypes;

  memory->create(tabindex,n+1,"bond:tabindex");
  memory->create(r0,n+1,"bond:r0");
  memory->create(setflag,n+1,"bond:setflag");
  for (int i = 1; i <= n; i++) setflag[i] = 0;
}


/* ----------------------------------------------------------------------
   global settings
------------------------------------------------------------------------- */

void BondTable::settings(int narg, char **arg)
{
  if (narg != 2) error->all(FLERR,"Illegal bond_style command");

  tabstyle = NONE;
  if (strcmp(arg[0],"linear") == 0) tabstyle = LINEAR;
  else if (strcmp(arg[0],"spline") == 0) tabstyle = SPLINE;
  else error->all(FLERR,"Unknown table style in bond style table");

  tablength = force->inumeric(FLERR,arg[1]);
  if (tablength < 2) error->all(FLERR,"Illegal number of bond table entries");

  // delete old tables, since cannot just change settings

  for (int m = 0; m < ntables; m++) free_table(&tables[m]);
  memory->sfree(tables);

  if (allocated) {
     memory->destroy(setflag);
     memory->destroy(tabindex);
  }
  allocated = 0;

  ntables = 0;
  tables = NULL;
}

/* ----------------------------------------------------------------------
   set coeffs for one or more type pairs
------------------------------------------------------------------------- */

void BondTable::coeff(int narg, char **arg)
{
  if (narg != 3) error->all(FLERR,"Illegal bond_coeff command");
  if (!allocated) allocate();

  int ilo,ihi;
  force->bounds(FLERR,arg[0],atom->nbondtypes,ilo,ihi);

  int me;
  MPI_Comm_rank(world,&me);
  tables = (Table *)
    memory->srealloc(tables,(ntables+1)*sizeof(Table),"bond:tables");
  Table *tb = &tables[ntables];
  null_table(tb);
  if (me == 0) read_table(tb,arg[1],arg[2]);
  bcast_table(tb);

  // error check on table parameters

  if (tb->ninput <= 1) error->one(FLERR,"Invalid bond table length");

  tb->lo = tb->rfile[0];
  tb->hi = tb->rfile[tb->ninput-1];
  if (tb->lo >= tb->hi) error->all(FLERR,"Bond table values are not increasing");

  // spline read-in and compute r,e,f vectors within table

  spline_table(tb);
  compute_table(tb);

  // store ptr to table in tabindex

  int count = 0;
  for (int i = ilo; i <= ihi; i++) {
    tabindex[i] = ntables;
    r0[i] = tb->r0;
    setflag[i] = 1;
    count++;
  }
  ntables++;

  if (count == 0) error->all(FLERR,"Illegal bond_coeff command");
}

/* ----------------------------------------------------------------------
   return an equilbrium bond length
   should not be used, since don't know minimum of tabulated function
------------------------------------------------------------------------- */

double BondTable::equilibrium_distance(int i)
{
  return r0[i];
}

/* ----------------------------------------------------------------------
   proc 0 writes to restart file
 ------------------------------------------------------------------------- */

void BondTable::write_restart(FILE *fp)
{
  write_restart_settings(fp);
}

/* ----------------------------------------------------------------------
    proc 0 reads from restart file, bcasts
 ------------------------------------------------------------------------- */

void BondTable::read_restart(FILE *fp)
{
  read_restart_settings(fp);
  allocate();
}

/* ----------------------------------------------------------------------
   proc 0 writes to restart file
 ------------------------------------------------------------------------- */

void BondTable::write_restart_settings(FILE *fp)
{
  fwrite(&tabstyle,sizeof(int),1,fp);
  fwrite(&tablength,sizeof(int),1,fp);
}

/* ----------------------------------------------------------------------
    proc 0 reads from restart file, bcasts
 ------------------------------------------------------------------------- */

void BondTable::read_restart_settings(FILE *fp)
{
  if (comm->me == 0) {
    fread(&tabstyle,sizeof(int),1,fp);
    fread(&tablength,sizeof(int),1,fp);
  }
  MPI_Bcast(&tabstyle,1,MPI_INT,0,world);
  MPI_Bcast(&tablength,1,MPI_INT,0,world);
}

/* ---------------------------------------------------------------------- */

double BondTable::single(int type, double rsq, int /*i*/, int /*j*/,
                         double &fforce)
{
  double r = sqrt(rsq);
  double u;
  double mdu;
  uf_lookup(type,r,u,mdu);
  fforce = mdu/r;
  return u;
}

/* ---------------------------------------------------------------------- */

void BondTable::null_table(Table *tb)
{
  tb->rfile = tb->efile = tb->ffile = NULL;
  tb->e2file = tb->f2file = NULL;
  tb->r = tb->e = tb->de = NULL;
  tb->f = tb->df = tb->e2 = tb->f2 = NULL;
}

/* ---------------------------------------------------------------------- */

void BondTable::free_table(Table *tb)
{
  memory->destroy(tb->rfile);
  memory->destroy(tb->efile);
  memory->destroy(tb->ffile);
  memory->destroy(tb->e2file);
  memory->destroy(tb->f2file);

  memory->destroy(tb->r);
  memory->destroy(tb->e);
  memory->destroy(tb->de);
  memory->destroy(tb->f);
  memory->destroy(tb->df);
  memory->destroy(tb->e2);
  memory->destroy(tb->f2);
}

/* ----------------------------------------------------------------------
   read table file, only called by proc 0
------------------------------------------------------------------------- */

void BondTable::read_table(Table *tb, char *file, char *keyword)
{
  char line[MAXLINE];
  double emin = BIGNUM;

  // open file

  FILE *fp = force->open_potential(file);
  if (fp == NULL) {
    char str[128];
    snprintf(str,128,"Cannot open file %s",file);
    error->one(FLERR,str);
  }

  // loop until section found with matching keyword

  while (1) {
    if (fgets(line,MAXLINE,fp) == NULL)
      error->one(FLERR,"Did not find keyword in table file");
    if (strspn(line," \t\n\r") == strlen(line)) continue;    // blank line
    if (line[0] == '#') continue;                          // comment
    char *word = strtok(line," \t\n\r");
    if (strcmp(word,keyword) == 0) break;            // matching keyword
    utils::sfgets(FLERR,line,MAXLINE,fp,file,error); // no match, skip section
    param_extract(tb,line);
    utils::sfgets(FLERR,line,MAXLINE,fp,file,error);
    for (int i = 0; i < tb->ninput; i++)
      utils::sfgets(FLERR,line,MAXLINE,fp,file,error);
  }

  // read args on 2nd line of section
  // allocate table arrays for file values

  utils::sfgets(FLERR,line,MAXLINE,fp,file,error);
  param_extract(tb,line);
  memory->create(tb->rfile,tb->ninput,"bond:rfile");
  memory->create(tb->efile,tb->ninput,"bond:efile");
  memory->create(tb->ffile,tb->ninput,"bond:ffile");

  // read r,e,f table values from file

  int cerror = 0;
  int r0idx = -1;

  utils::sfgets(FLERR,line,MAXLINE,fp,file,error);
  for (int i = 0; i < tb->ninput; i++) {
    if (NULL == fgets(line,MAXLINE,fp))
      error->one(FLERR,"Premature end of file in bond table");
    if (3 != sscanf(line,"%*d %lg %lg %lg",
                    &tb->rfile[i],&tb->efile[i],&tb->ffile[i])) ++cerror;
    if (tb->efile[i] < emin) {
      emin = tb->efile[i];
      r0idx = i;
    }
  }
  fclose(fp);

  // infer r0 from minimum of potential, if not given explicitly

  if ((tb->r0 == 0.0) && (r0idx >= 0)) tb->r0 = tb->rfile[r0idx];

  // warn if force != dE/dr at any point that is not an inflection point
  // check via secant approximation to dE/dr
  // skip two end points since do not have surrounding secants
  // inflection point is where curvature changes sign

  double r,e,f,rprev,rnext,eprev,enext,fleft,fright;

  int ferror = 0;
  for (int i = 1; i < tb->ninput-1; i++) {
    r = tb->rfile[i];
    rprev = tb->rfile[i-1];
    rnext = tb->rfile[i+1];
    e = tb->efile[i];
    eprev = tb->efile[i-1];
    enext = tb->efile[i+1];
    f = tb->ffile[i];
    fleft = - (e-eprev) / (r-rprev);
    fright = - (enext-e) / (rnext-r);
    if (f < fleft && f < fright) ferror++;
    if (f > fleft && f > fright) ferror++;
  }

  if (ferror) {
    char str[128];
    sprintf(str,"%d of %d force values in table are inconsistent with -dE/dr.\n"
            "  Should only be flagged at inflection points",ferror,tb->ninput);
    error->warning(FLERR,str);
  }

  // warn if data was read incompletely, e.g. columns were missing

  if (cerror) {
    char str[128];
    sprintf(str,"%d of %d lines in table were incomplete or could not be"
            " parsed completely",cerror,tb->ninput);
    error->warning(FLERR,str);
  }
}

/* ----------------------------------------------------------------------
   build spline representation of e,f over entire range of read-in table
   this function sets these values in e2file,f2file
------------------------------------------------------------------------- */

void BondTable::spline_table(Table *tb)
{
  memory->create(tb->e2file,tb->ninput,"bond:e2file");
  memory->create(tb->f2file,tb->ninput,"bond:f2file");

  double ep0 = - tb->ffile[0];
  double epn = - tb->ffile[tb->ninput-1];
  spline(tb->rfile,tb->efile,tb->ninput,ep0,epn,tb->e2file);

  if (tb->fpflag == 0) {
    tb->fplo = (tb->ffile[1] - tb->ffile[0]) / (tb->rfile[1] - tb->rfile[0]);
    tb->fphi = (tb->ffile[tb->ninput-1] - tb->ffile[tb->ninput-2]) /
      (tb->rfile[tb->ninput-1] - tb->rfile[tb->ninput-2]);
  }

  double fp0 = tb->fplo;
  double fpn = tb->fphi;
  spline(tb->rfile,tb->ffile,tb->ninput,fp0,fpn,tb->f2file);

}

/* ----------------------------------------------------------------------
   compute r,e,f vectors from splined values
------------------------------------------------------------------------- */

void BondTable::compute_table(Table *tb)
{
  // delta = table spacing for N-1 bins
  int tlm1 = tablength-1;

  tb->delta = (tb->hi - tb->lo)/ tlm1;
  tb->invdelta = 1.0/tb->delta;
  tb->deltasq6 = tb->delta*tb->delta / 6.0;

  // N-1 evenly spaced bins in r from min to max
  // r,e,f = value at lower edge of bin
  // de,df values = delta values of e,f
  // r,e,f are N in length so de,df arrays can compute difference

  memory->create(tb->r,tablength,"bond:r");
  memory->create(tb->e,tablength,"bond:e");
  memory->create(tb->de,tlm1,"bond:de");
  memory->create(tb->f,tablength,"bond:f");
  memory->create(tb->df,tlm1,"bond:df");
  memory->create(tb->e2,tablength,"bond:e2");
  memory->create(tb->f2,tablength,"bond:f2");

  double a;
  for (int i = 0; i < tablength; i++) {
    a = tb->lo + i*tb->delta;
    tb->r[i] = a;
    tb->e[i] = splint(tb->rfile,tb->efile,tb->e2file,tb->ninput,a);
    tb->f[i] = splint(tb->rfile,tb->ffile,tb->f2file,tb->ninput,a);
  }

  for (int i = 0; i < tlm1; i++) {
    tb->de[i] = tb->e[i+1] - tb->e[i];
    tb->df[i] = tb->f[i+1] - tb->f[i];
  }

  double ep0 = - tb->f[0];
  double epn = - tb->f[tlm1];
  spline(tb->r,tb->e,tablength,ep0,epn,tb->e2);
  spline(tb->r,tb->f,tablength,tb->fplo,tb->fphi,tb->f2);
}

/* ----------------------------------------------------------------------
   extract attributes from parameter line in table section
   format of line: N value FP fplo fphi EQ r0
   N is required, other params are optional
------------------------------------------------------------------------- */

void BondTable::param_extract(Table *tb, char *line)
{
  tb->ninput = 0;
  tb->fpflag = 0;
  tb->r0 = 0.0;

  char *word = strtok(line," \t\n\r\f");
  while (word) {
    if (strcmp(word,"N") == 0) {
      word = strtok(NULL," \t\n\r\f");
      tb->ninput = atoi(word);
    } else if (strcmp(word,"FP") == 0) {
      tb->fpflag = 1;
      word = strtok(NULL," \t\n\r\f");
      tb->fplo = atof(word);
      word = strtok(NULL," \t\n\r\f");
      tb->fphi = atof(word);
    } else if (strcmp(word,"EQ") == 0) {
      word = strtok(NULL," \t\n\r\f");
      tb->r0 = atof(word);
    } else {
      error->one(FLERR,"Invalid keyword in bond table parameters");
    }
    word = strtok(NULL," \t\n\r\f");
  }

  if (tb->ninput == 0) error->one(FLERR,"Bond table parameters did not set N");
}

/* ----------------------------------------------------------------------
   broadcast read-in table info from proc 0 to other procs
   this function communicates these values in Table:
     ninput,rfile,efile,ffile,fpflag,fplo,fphi,r0
------------------------------------------------------------------------- */

void BondTable::bcast_table(Table *tb)
{
  MPI_Bcast(&tb->ninput,1,MPI_INT,0,world);
  MPI_Bcast(&tb->r0,1,MPI_DOUBLE,0,world);

  int me;
  MPI_Comm_rank(world,&me);
  if (me > 0) {
    memory->create(tb->rfile,tb->ninput,"angle:rfile");
    memory->create(tb->efile,tb->ninput,"angle:efile");
    memory->create(tb->ffile,tb->ninput,"angle:ffile");
  }

  MPI_Bcast(tb->rfile,tb->ninput,MPI_DOUBLE,0,world);
  MPI_Bcast(tb->efile,tb->ninput,MPI_DOUBLE,0,world);
  MPI_Bcast(tb->ffile,tb->ninput,MPI_DOUBLE,0,world);

  MPI_Bcast(&tb->fpflag,1,MPI_INT,0,world);
  if (tb->fpflag) {
    MPI_Bcast(&tb->fplo,1,MPI_DOUBLE,0,world);
    MPI_Bcast(&tb->fphi,1,MPI_DOUBLE,0,world);
  }
}

/* ----------------------------------------------------------------------
   spline and splint routines modified from Numerical Recipes
------------------------------------------------------------------------- */

void BondTable::spline(double *x, double *y, int n,
                       double yp1, double ypn, double *y2)
{
  int i,k;
  double p,qn,sig,un;
  double *u = new double[n];

  if (yp1 > 0.99e30) y2[0] = u[0] = 0.0;
  else {
    y2[0] = -0.5;
    u[0] = (3.0/(x[1]-x[0])) * ((y[1]-y[0]) / (x[1]-x[0]) - yp1);
  }
  for (i = 1; i < n-1; i++) {
    sig = (x[i]-x[i-1]) / (x[i+1]-x[i-1]);
    p = sig*y2[i-1] + 2.0;
    y2[i] = (sig-1.0) / p;
    u[i] = (y[i+1]-y[i]) / (x[i+1]-x[i]) - (y[i]-y[i-1]) / (x[i]-x[i-1]);
    u[i] = (6.0*u[i] / (x[i+1]-x[i-1]) - sig*u[i-1]) / p;
  }
  if (ypn > 0.99e30) qn = un = 0.0;
  else {
    qn = 0.5;
    un = (3.0/(x[n-1]-x[n-2])) * (ypn - (y[n-1]-y[n-2]) / (x[n-1]-x[n-2]));
  }
  y2[n-1] = (un-qn*u[n-2]) / (qn*y2[n-2] + 1.0);
  for (k = n-2; k >= 0; k--) y2[k] = y2[k]*y2[k+1] + u[k];

  delete [] u;
}

/* ---------------------------------------------------------------------- */

double BondTable::splint(double *xa, double *ya, double *y2a, int n, double x)
{
  int klo,khi,k;
  double h,b,a,y;

  klo = 0;
  khi = n-1;
  while (khi-klo > 1) {
    k = (khi+klo) >> 1;
    if (xa[k] > x) khi = k;
    else klo = k;
  }
  h = xa[khi]-xa[klo];
  a = (xa[khi]-x) / h;
  b = (x-xa[klo]) / h;
  y = a*ya[klo] + b*ya[khi] +
    ((a*a*a-a)*y2a[klo] + (b*b*b-b)*y2a[khi]) * (h*h)/6.0;
  return y;
}

/* ----------------------------------------------------------------------
   calculate potential u and force f at distance x
   insure x is between bond min/max, exit with error if not
------------------------------------------------------------------------- */

void BondTable::uf_lookup(int type, double x, double &u, double &f)
{
  if (!std::isfinite(x)) {
    error->one(FLERR,"Illegal bond in bond style table");
  }

  double fraction,a,b;
  char estr[128];
  const Table *tb = &tables[tabindex[type]];
  const int itable = static_cast<int> ((x - tb->lo) * tb->invdelta);
  if (itable < 0) {
    sprintf(estr,"Bond length < table inner cutoff: "
            "type %d length %g",type,x);
    error->one(FLERR,estr);
  } else if (itable >= tablength) {
    sprintf(estr,"Bond length > table outer cutoff: "
            "type %d length %g",type,x);
    error->one(FLERR,estr);
  }

  if (tabstyle == LINEAR) {
    fraction = (x - tb->r[itable]) * tb->invdelta;
    u = tb->e[itable] + fraction*tb->de[itable];
    f = tb->f[itable] + fraction*tb->df[itable];
  } else if (tabstyle == SPLINE) {
    fraction = (x - tb->r[itable]) * tb->invdelta;

    b = (x - tb->r[itable]) * tb->invdelta;
    a = 1.0 - b;
    u = a * tb->e[itable] + b * tb->e[itable+1] +
      ((a*a*a-a)*tb->e2[itable] + (b*b*b-b)*tb->e2[itable+1]) *
      tb->deltasq6;
    f = a * tb->f[itable] + b * tb->f[itable+1] +
      ((a*a*a-a)*tb->f2[itable] + (b*b*b-b)*tb->f2[itable+1]) *
      tb->deltasq6;
  }
}

/* ----------------------------------------------------------------------
   calculate potential u at distance x
   insure x is between bond min/max
------------------------------------------------------------------------- */

void BondTable::u_lookup(int type, double x, double &u)
{
  if (!std::isfinite(x)) {
    error->one(FLERR,"Illegal bond in bond style table");
  }

  double fraction,a,b;
  char estr[128];
  const Table *tb = &tables[tabindex[type]];
  const int itable = static_cast<int> ((x - tb->lo) * tb->invdelta);
  if (itable < 0) {
    sprintf(estr,"Bond length < table inner cutoff: "
            "type %d length %g",type,x);
    error->one(FLERR,estr);
  } else if (itable >= tablength) {
    sprintf(estr,"Bond length > table outer cutoff: "
            "type %d length %g",type,x);
    error->one(FLERR,estr);
  }

  if (tabstyle == LINEAR) {
    fraction = (x - tb->r[itable]) * tb->invdelta;
    u = tb->e[itable] + fraction*tb->de[itable];
  } else if (tabstyle == SPLINE) {
    fraction = (x - tb->r[itable]) * tb->invdelta;

    b = (x - tb->r[itable]) * tb->invdelta;
    a = 1.0 - b;
    u = a * tb->e[itable] + b * tb->e[itable+1] +
      ((a*a*a-a)*tb->e2[itable] + (b*b*b-b)*tb->e2[itable+1]) *
      tb->deltasq6;
  }
}