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
/* ----------------------------------------------------------------------
   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 authors: Stephen Bond (SNL) and
                         Andrew Baczewski (Michigan State/SNL)
------------------------------------------------------------------------- */

#include "fix_gld.h"
#include <mpi.h>
#include <cmath>
#include <cstring>
#include "atom.h"
#include "force.h"
#include "update.h"
#include "respa.h"
#include "comm.h"
#include "random_mars.h"
#include "memory.h"
#include "error.h"
#include "group.h"

#define GLD_UNIFORM_DISTRO

using namespace LAMMPS_NS;
using namespace FixConst;

/* ----------------------------------------------------------------------
   Parses parameters passed to the method, allocates some memory
------------------------------------------------------------------------- */

FixGLD::FixGLD(LAMMPS *lmp, int narg, char **arg) :
  Fix(lmp, narg, arg),
  step_respa(NULL), prony_c(NULL), prony_tau(NULL), s_gld(NULL), random(NULL)
{
  int narg_min = 8;
  // Check to make sure we have the minimal number of inputs
  if (narg < narg_min) error->all(FLERR,"Illegal fix gld command");

  time_integrate = 1;
  restart_peratom = 1;

  // Parse the first set of required input arguments
  // 0 = Fix ID           (e.g., 1)
  // 1 = Group ID         (e.g., all)
  // 2 = gld              (name of this fix)
  // 3 = t_start          (Starting target temperature)
  t_start = force->numeric(FLERR,arg[3]);

  // 4 = t_stop           (Stopping target temperature)
  t_stop = force->numeric(FLERR,arg[4]);

  // 5 = prony_terms      (number of terms in Prony series)
  prony_terms = force->inumeric(FLERR,arg[5]);

  // 6 = seed             (random seed)
  int seed    = force->inumeric(FLERR,arg[6]);

  // 7 = series type
  if(strcmp(arg[7],"pprony") == 0) {
     series_type = 1;   // series type 1 is 'positive Prony series'
  } else {
     error->all(FLERR,"Fix gld series type must be pprony for now");
  }

  // Error checking for the first set of required input arguments
  if (seed <= 0) error->all(FLERR,"Illegal fix gld command");
  if (prony_terms <= 0)
    error->all(FLERR,"Fix gld prony terms must be > 0");
  if (t_start < 0)
    error->all(FLERR,"Fix gld start temperature must be >= 0");
  if (t_stop < 0)
    error->all(FLERR,"Fix gld stop temperature must be >= 0");
  if (narg - narg_min < 2*(prony_terms) )
    error->all(FLERR,"Fix gld needs more prony series coefficients");

  // allocate memory for Prony series force coefficients
  memory->create(prony_c, prony_terms, "gld:prony_c");
  // allocate memory for Prony series timescale coefficients
  memory->create(prony_tau, prony_terms, "gld:prony_tau");
  // allocate memory for Prony series extended variables
  s_gld = NULL;
  grow_arrays(atom->nmax);
  // add callbacks to enable restarts
  atom->add_callback(0);
  atom->add_callback(1);

  // read in the Prony series coefficients
  int iarg = narg_min;
  int icoeff = 0;
  while (iarg < narg && icoeff < prony_terms) {
    double pc = force->numeric(FLERR,arg[iarg]);
    double ptau = force->numeric(FLERR,arg[iarg+1]);

    if (pc < 0)
      error->all(FLERR,"Fix gld c coefficients must be >= 0");
    if (ptau <= 0)
      error->all(FLERR,"Fix gld tau coefficients must be > 0");

    // All atom types to have the same Prony series
    prony_c[icoeff] = pc;
    prony_tau[icoeff] = ptau;

    icoeff += 1;
    iarg += 2;
  }

  // initialize Marsaglia RNG with processor-unique seed
  random = new RanMars(lmp,seed + comm->me);

  // initialize the extended variables
  init_s_gld();

  // optional arguments
  freezeflag = 0;
  zeroflag = 0;

  while (iarg < narg) {
    if (strcmp(arg[iarg],"zero") == 0) {
      if (iarg+2 > narg) {
        error->all(FLERR, "Illegal fix gld command");
      }
      if (strcmp(arg[iarg+1],"no") == 0) {
      } else if (strcmp(arg[iarg+1],"yes") == 0) {
        zeroflag = 1;
      } else {
        error->all(FLERR,"Illegal fix gld command");
      }
      iarg += 2;
    }
    else if (strcmp(arg[iarg],"frozen") == 0) {
       if (iarg+2 > narg) {
          error->all(FLERR, "Illegal fix gld command");
       }
       if (strcmp(arg[iarg+1],"no") == 0) {
       } else if (strcmp(arg[iarg+1],"yes") == 0) {
         freezeflag = 1;
         for (int i = 0; i < atom->nlocal; i++) {
           if (atom->mask[i] & groupbit) {
             for (int k = 0; k < 3*prony_terms; k=k+3)
             {
               s_gld[i][k] = 0.0;
               s_gld[i][k+1] = 0.0;
               s_gld[i][k+2] = 0.0;
             }
           }
         }
       } else {
          error->all(FLERR, "Illegal fix gld command");
       }
       iarg += 2;
    }
    else error->all(FLERR,"Illegal fix gld command");
  }

  // Initialize the target temperature
  t_target = t_start;
}

/* ----------------------------------------------------------------------
   Destroys memory allocated by the method
------------------------------------------------------------------------- */

FixGLD::~FixGLD()
{
  delete random;
  memory->destroy(prony_c);
  memory->destroy(prony_tau);
  memory->destroy(s_gld);

  // remove callbacks to fix, so atom class stops calling it
  atom->delete_callback(id,0);
  atom->delete_callback(id,1);
}

/* ----------------------------------------------------------------------
   Specifies when the fix is called during the timestep
------------------------------------------------------------------------- */

int FixGLD::setmask()
{
  int mask = 0;
  mask |= INITIAL_INTEGRATE;
  mask |= FINAL_INTEGRATE;
  mask |= INITIAL_INTEGRATE_RESPA;
  mask |= FINAL_INTEGRATE_RESPA;
  return mask;
}

/* ----------------------------------------------------------------------
   Initialize the method parameters before a run
------------------------------------------------------------------------- */

void FixGLD::init()
{
  dtv = update->dt;
  dtf = 0.5 * update->dt * force->ftm2v;

  if (strstr(update->integrate_style,"respa"))
    step_respa = ((Respa *) update->integrate)->step;
}

/* ----------------------------------------------------------------------
   First half of a timestep (V^{n} -> V^{n+1/2}; X^{n} -> X^{n+1})
------------------------------------------------------------------------- */

void FixGLD::initial_integrate(int /*vflag*/)
{
  double dtfm;
  double ftm2v = force->ftm2v;

  double fran[3], fsum[3], fsumall[3];
  bigint count;

  int icoeff;

  // update v and x of atoms in group
  double **x = atom->x;
  double **v = atom->v;
  double **f = atom->f;
  double *rmass = atom->rmass;
  double *mass = atom->mass;
  int *type = atom->type;
  int *mask = atom->mask;
  int nlocal = atom->nlocal;
  if (igroup == atom->firstgroup) nlocal = atom->nfirst;

  // set kT to the temperature in mvvv units
  double kT = (force->boltz)*t_target/(force->mvv2e);

  // zero an accumulator for the total random force
  fsum[0] = fsum[1] = fsum[2] = 0.0;

  if (rmass) {
    for (int i = 0; i < nlocal; i++)
      if (mask[i] & groupbit) {
        dtfm = dtf / rmass[i];
        // Advance V by dt/2
        v[i][0] += dtfm * f[i][0];
        v[i][1] += dtfm * f[i][1];
        v[i][2] += dtfm * f[i][2];
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          v[i][0] += dtfm * s_gld[i][k];
          v[i][1] += dtfm * s_gld[i][k+1];
          v[i][2] += dtfm * s_gld[i][k+2];
        }

        // Advance X by dt
        x[i][0] += dtv * v[i][0];
        x[i][1] += dtv * v[i][1];
        x[i][2] += dtv * v[i][2];

        // Advance S by dt
        icoeff = 0;
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          double theta = exp(-dtv/prony_tau[icoeff]);
          double ck = prony_c[icoeff];
          double vmult = (theta-1.)*ck/ftm2v;
          double rmult = sqrt(2.0*kT*ck/dtv)*(1.-theta)/ftm2v;

          // random force
#ifdef GLD_GAUSSIAN_DISTRO
          fran[0] = rmult*random->gaussian();
          fran[1] = rmult*random->gaussian();
          fran[2] = rmult*random->gaussian();
#endif

#ifdef GLD_UNIFORM_DISTRO
          rmult *= sqrt(12.0); // correct variance of uniform distribution
          fran[0] = rmult*(random->uniform() - 0.5);
          fran[1] = rmult*(random->uniform() - 0.5);
          fran[2] = rmult*(random->uniform() - 0.5);
#endif

          // sum of random forces
          fsum[0] += fran[0];
          fsum[1] += fran[1];
          fsum[2] += fran[2];

          s_gld[i][k]   *= theta;
          s_gld[i][k+1] *= theta;
          s_gld[i][k+2] *= theta;
          s_gld[i][k]   += vmult*v[i][0];
          s_gld[i][k+1] += vmult*v[i][1];
          s_gld[i][k+2] += vmult*v[i][2];
          s_gld[i][k]   += fran[0];
          s_gld[i][k+1] += fran[1];
          s_gld[i][k+2] += fran[2];

          icoeff += 1;
        }
      }

  } else {
    for (int i = 0; i < nlocal; i++)
      if (mask[i] & groupbit) {
        dtfm = dtf / mass[type[i]];
        // Advance V by dt/2
        v[i][0] += dtfm * f[i][0];
        v[i][1] += dtfm * f[i][1];
        v[i][2] += dtfm * f[i][2];
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          v[i][0] += dtfm * s_gld[i][k];
          v[i][1] += dtfm * s_gld[i][k+1];
          v[i][2] += dtfm * s_gld[i][k+2];
        }

        // Advance X by dt
        x[i][0] += dtv * v[i][0];
        x[i][1] += dtv * v[i][1];
        x[i][2] += dtv * v[i][2];

        // Advance S by dt
        icoeff = 0;
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          double theta = exp(-dtv/prony_tau[icoeff]);
          double ck = prony_c[icoeff];
          double vmult = (theta-1.)*ck/ftm2v;
          double rmult = sqrt(2.0*kT*ck/dtv)*(1.-theta)/ftm2v;

          // random force
#ifdef GLD_GAUSSIAN_DISTRO
          fran[0] = rmult*random->gaussian();
          fran[1] = rmult*random->gaussian();
          fran[2] = rmult*random->gaussian();
#endif

#ifdef GLD_UNIFORM_DISTRO
          rmult *= sqrt(12.0); // correct variance of uniform distribution
          fran[0] = rmult*(random->uniform() - 0.5);
          fran[1] = rmult*(random->uniform() - 0.5);
          fran[2] = rmult*(random->uniform() - 0.5);
#endif

          // sum of random forces
          fsum[0] += fran[0];
          fsum[1] += fran[1];
          fsum[2] += fran[2];

          s_gld[i][k]   *= theta;
          s_gld[i][k+1] *= theta;
          s_gld[i][k+2] *= theta;
          s_gld[i][k]   += vmult*v[i][0];
          s_gld[i][k+1] += vmult*v[i][1];
          s_gld[i][k+2] += vmult*v[i][2];
          s_gld[i][k]   += fran[0];
          s_gld[i][k+1] += fran[1];
          s_gld[i][k+2] += fran[2];

          icoeff += 1;

        }
      }
  }

  // correct the random force, if zeroflag is set
  if (zeroflag) {
    count = group->count(igroup);
    if (count == 0) error->all(FLERR,"Cannot zero gld force for zero atoms");

    MPI_Allreduce(fsum,fsumall,3,MPI_DOUBLE,MPI_SUM,world);
    fsumall[0] /= (count*prony_terms);
    fsumall[1] /= (count*prony_terms);
    fsumall[2] /= (count*prony_terms);
    for (int i = 0; i < nlocal; i++) {
      if (mask[i] & groupbit) {
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          s_gld[i][k]   -= fsumall[0];
          s_gld[i][k+1] -= fsumall[1];
          s_gld[i][k+2] -= fsumall[2];
        }
      }
    }
  }

}

/* ----------------------------------------------------------------------
   Second half of a timestep (V^{n+1/2} -> V^{n+1})
------------------------------------------------------------------------- */

void FixGLD::final_integrate()
{
  double dtfm;

  // update v of atoms in group

  double **v = atom->v;
  double **f = atom->f;
  double *rmass = atom->rmass;
  double *mass = atom->mass;
  int *type = atom->type;
  int *mask = atom->mask;
  int nlocal = atom->nlocal;
  if (igroup == atom->firstgroup) nlocal = atom->nfirst;

  if (rmass) {
    for (int i = 0; i < nlocal; i++)
      if (mask[i] & groupbit) {
        dtfm = dtf / rmass[i];
        v[i][0] += dtfm * f[i][0];
        v[i][1] += dtfm * f[i][1];
        v[i][2] += dtfm * f[i][2];
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          v[i][0] += dtfm * s_gld[i][k];
          v[i][1] += dtfm * s_gld[i][k+1];
          v[i][2] += dtfm * s_gld[i][k+2];
        }
      }

  } else {
    for (int i = 0; i < nlocal; i++)
      if (mask[i] & groupbit) {
        dtfm = dtf / mass[type[i]];
        v[i][0] += dtfm * f[i][0];
        v[i][1] += dtfm * f[i][1];
        v[i][2] += dtfm * f[i][2];
        for (int k = 0; k < 3*prony_terms; k=k+3) {
          v[i][0] += dtfm * s_gld[i][k];
          v[i][1] += dtfm * s_gld[i][k+1];
          v[i][2] += dtfm * s_gld[i][k+2];
        }
      }
  }

  // Change the temperature for the next step
  double delta = update->ntimestep - update->beginstep;
  delta /= update->endstep - update->beginstep;
  t_target = t_start + delta * (t_stop - t_start);
}

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

void FixGLD::initial_integrate_respa(int vflag, int ilevel, int /*iloop*/)
{
  dtv = step_respa[ilevel];
  dtf = 0.5 * step_respa[ilevel] * (force->ftm2v);

  // innermost level - GLD update of v and x
  // all other levels - GLD update of v

  if (ilevel == 0) initial_integrate(vflag);
  else final_integrate();
}

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

void FixGLD::final_integrate_respa(int ilevel, int /*iloop*/)
{
  dtf = 0.5 * step_respa[ilevel] * (force->ftm2v);
  final_integrate();
}

/* ----------------------------------------------------------------------
   Called when a change to the target temperature is requested mid-run
------------------------------------------------------------------------- */

void FixGLD::reset_target(double t_new)
{
  t_target = t_start = t_stop = t_new;
}

/* ----------------------------------------------------------------------
   Called when a change to the timestep is requested mid-run
------------------------------------------------------------------------- */

void FixGLD::reset_dt()
{
  // set the time integration constants
  dtv = update->dt;
  dtf = 0.5 * update->dt * (force->ftm2v);
}

/* ----------------------------------------------------------------------
   memory usage of local atom-based arrays
------------------------------------------------------------------------- */

double FixGLD::memory_usage()
{
  double bytes = atom->nmax*3*prony_terms*sizeof(double);
  return bytes;
}

/* ----------------------------------------------------------------------
   allocate local atom-based arrays
------------------------------------------------------------------------- */

void FixGLD::grow_arrays(int nmax)
{
  memory->grow(s_gld, nmax, 3*prony_terms,"gld:s_gld");
}

/* ----------------------------------------------------------------------
   copy values within local atom-based arrays
------------------------------------------------------------------------- */

void FixGLD::copy_arrays(int i, int j, int /*delflag*/)
{
  for (int k = 0; k < 3*prony_terms; k++) {
    s_gld[j][k] = s_gld[i][k];
  }
}

/* ----------------------------------------------------------------------
   Pack extended variables assoc. w/ atom i into buffer for exchange
   with another processor
------------------------------------------------------------------------- */

int FixGLD::pack_exchange(int i, double *buf)
{
  int m = 0;
  for (int k = 0; k < 3*prony_terms; k++) {
    buf[m++] = s_gld[i][k];
  }
  return m;
}

/* ----------------------------------------------------------------------
   Unpack extended variables from exchange with another processor
------------------------------------------------------------------------- */

int FixGLD::unpack_exchange(int nlocal, double *buf)
{
  int m = 0;
  for (int k = 0; k < 3*prony_terms; k++) {
    s_gld[nlocal][k] = buf[m++];
  }
  return m;
}


/* ----------------------------------------------------------------------
   Pack extended variables assoc. w/ atom i into buffer for
   writing to a restart file
------------------------------------------------------------------------- */

int FixGLD::pack_restart(int i, double *buf)
{
  int m = 0;
  buf[m++] = 3*prony_terms + 1;
  for (int k = 0; k < 3*prony_terms; k=k+3)
  {
    buf[m++] = s_gld[i][k];
    buf[m++] = s_gld[i][k+1];
    buf[m++] = s_gld[i][k+2];
  }
  return m;
}

/* ----------------------------------------------------------------------
   Unpack extended variables to restart the fix from a restart file
------------------------------------------------------------------------- */

void FixGLD::unpack_restart(int nlocal, int nth)
{
  double **extra = atom->extra;

  // skip to the nth set of extended variables

  int m = 0;
  for (int i = 0; i< nth; i++) m += static_cast<int> (extra[nlocal][m]);
  m++;

  for (int k = 0; k < 3*prony_terms; k=k+3)
  {
    s_gld[nlocal][k] = extra[nlocal][m++];
    s_gld[nlocal][k+1] = extra[nlocal][m++];
    s_gld[nlocal][k+2] = extra[nlocal][m++];
  }
}

/* ----------------------------------------------------------------------
   Returns the number of items in atomic restart data associated with
   local atom nlocal.  Used in determining the total extra data stored by
   fixes on a given processor.
------------------------------------------------------------------------- */

int FixGLD::size_restart(int /*nlocal*/)
{
  return 3*prony_terms+1;
}

/* ----------------------------------------------------------------------
   Returns the maximum number of items in atomic restart data
   Called in Modify::restart for peratom restart.
------------------------------------------------------------------------- */

int FixGLD::maxsize_restart()
{
  return 3*prony_terms+1;
}

/* ----------------------------------------------------------------------
   Initializes the extended variables to equilibrium distribution
   at t_start.
------------------------------------------------------------------------- */

void FixGLD::init_s_gld()
{
  int icoeff;
  double eq_sdev=0.0;

  // set kT to the temperature in mvvv units
  double kT = (force->boltz)*t_start/(force->mvv2e);

#ifdef GLD_GAUSSIAN_DISTRO
  double scale = sqrt(kT)/(force->ftm2v);
#endif

#ifdef GLD_UNIFORM_DISTRO
  double scale = sqrt(12.0*kT)/(force->ftm2v);
#endif

  for (int i = 0; i < atom->nlocal; i++) {
    if (atom->mask[i] & groupbit) {
      icoeff = 0;
      for (int k = 0; k < 3*prony_terms; k=k+3) {
        eq_sdev = scale*sqrt(prony_c[icoeff]/prony_tau[icoeff]);
#ifdef GLD_GAUSSIAN_DISTRO
        s_gld[i][k] = eq_sdev*random->gaussian();
        s_gld[i][k+1] = eq_sdev*random->gaussian();
        s_gld[i][k+2] = eq_sdev*random->gaussian();
#endif

#ifdef GLD_UNIFORM_DISTRO
        s_gld[i][k] = eq_sdev*(random->uniform()-0.5);
        s_gld[i][k+1] = eq_sdev*(random->uniform()-0.5);
        s_gld[i][k+2] = eq_sdev*(random->uniform()-0.5);
#endif
        icoeff += 1;
      }
    }
  }

  return;
}