File: mc.py

package info (click to toggle)
lammps 20250204%2Bdfsg.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 474,368 kB
  • sloc: cpp: 1,060,070; python: 27,785; ansic: 8,956; f90: 7,254; sh: 6,044; perl: 4,171; fortran: 2,442; xml: 1,714; makefile: 1,352; objc: 238; lisp: 188; yacc: 58; csh: 16; awk: 14; tcl: 6; javascript: 2
file content (109 lines) | stat: -rwxr-xr-x 2,724 bytes parent folder | download
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
#!/usr/bin/env python -i
# preceding line should have path for Python on your machine

# mc.py
# Purpose: mimic operation of example/MC/in.mc via Python
# Syntax:  mc.py in.mc
#          in.mc = LAMMPS input script

from __future__ import print_function
import sys,random,math

# set these parameters
# make sure neigh skin (in in.mc) > 2*deltamove

nloop = 3000
deltaperturb = 0.2
deltamove = 0.1
kT = 0.05
random.seed(27848)

# parse command line

argv = sys.argv
if len(argv) != 2:
  print("Syntax: mc.py in.mc")
  sys.exit()

infile = sys.argv[1]

from lammps import lammps, LAMMPS_INT, LMP_STYLE_GLOBAL, LMP_VAR_EQUAL
lmp = lammps()

# run infile one line at a time
# just sets up MC problem

lines = open(infile,'r').readlines()
for line in lines: lmp.command(line)
lmp.command("variable e equal pe")

# run 0 to get energy of perfect lattice
# emin = minimum energy

lmp.command("run 0")

natoms = lmp.extract_global("natoms")
emin = lmp.extract_compute("thermo_pe",LMP_STYLE_GLOBAL,LAMMPS_INT) / natoms
lmp.command("variable emin equal $e")

# disorder the system
# estart = initial energy

x = lmp.extract_atom("x")

for i in range(natoms):
  x[i][0] += deltaperturb * (2*random.random()-1)
  x[i][1] += deltaperturb * (2*random.random()-1)

lmp.command("variable elast equal $e")
lmp.command("thermo_style custom step v_emin v_elast pe")
lmp.command("run 0")
x = lmp.extract_atom("x")
lmp.command("variable elast equal $e")

estart = lmp.extract_compute("thermo_pe", LMP_STYLE_GLOBAL, LAMMPS_INT) / natoms

# loop over Monte Carlo moves
# extract x after every run, in case reneighboring changed ptr in LAMMPS

elast = estart
naccept = 0

for i in range(nloop):
  iatom = random.randrange(0,natoms)
  x0 = x[iatom][0]
  y0 = x[iatom][1]

  x[iatom][0] += deltamove * (2*random.random()-1)
  x[iatom][1] += deltamove * (2*random.random()-1)

  lmp.command("run 1 pre no post no")
  x = lmp.extract_atom("x")
  e = lmp.extract_compute("thermo_pe", LMP_STYLE_GLOBAL, LAMMPS_INT) / natoms

  if e <= elast:
    elast = e
    lmp.command("variable elast equal $e")
    naccept += 1
  elif random.random() <= math.exp(natoms*(elast-e)/kT):
    elast = e
    lmp.command("variable elast equal $e")
    naccept += 1
  else:
    x[iatom][0] = x0
    x[iatom][1] = y0

# final energy and stats

lmp.command("variable nbuild equal nbuild")
nbuild = lmp.extract_variable("nbuild", None, LMP_VAR_EQUAL)

lmp.command("run 0")
estop = lmp.extract_compute("thermo_pe", LMP_STYLE_GLOBAL, LAMMPS_INT) / natoms

print("MC stats:")
print("  starting energy =",estart)
print("  final energy =",estop)
print("  minimum energy of perfect lattice =",emin)
print("  accepted MC moves =",naccept)
print("  neighbor list rebuilds =",nbuild)