-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathSurface.py
145 lines (126 loc) · 6.05 KB
/
Surface.py
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
from Property import Property
from refine import make_refine
import reproduce
import dpgen.auto_test.lib.vasp as vasp
import dpgen.auto_test.lib.lammps as lammps
from pymatgen.core.structure import Structure
from pymatgen.core.surface import generate_all_slabs
import numpy as np
import os,json
class Surface(Property):
def __init__(self,
parameter):
self.parameter = parameter
self.min_slab_size = parameter['min_slab_size']
self.min_vacuum_size = parameter['min_vacuum_size']
self.pert_xz = parameter['pert_xz']
default_max_miller = 2
self.miller = parameter.get('max_miller', default_max_miller)
self.static = parameter.get('static-opt', False)
self.relax = parameter.get('change_box', False)
self.reprod = parameter.get('reprod-opt', False)
def make_confs(self,
path_to_work,
path_to_equi,
refine=False):
path_to_work = os.path.abspath(path_to_work)
path_to_equi = os.path.abspath(path_to_equi)
task_list = []
cwd = os.getcwd()
equi_contcar = os.path.join(path_to_equi, 'CONTCAR')
ptypes = vasp.get_poscar_types(equi_contcar)
if not os.path.exists(equi_contcar):
raise RuntimeError("please do relaxation first")
# gen structure
ss = Structure.from_file(equi_contcar)
# gen slabs
all_slabs = generate_all_slabs(ss, self.max_miller, self.min_slab_size, self.min_vacuum_size)
if refine:
task_list = make_refine(self.parameter['init_from_suffix'],
self.parameter['output_suffix'],
path_to_work,
len(all_slabs))
# record miller
for ii in range(len(task_list)):
os.chdir(task_list[ii])
np.savetxt('miller.out', all_slabs[ii].miller_index, fmt='%d')
os.chdir(cwd)
if self.reprod:
if 'vasp_path' not in self.parameter:
raise RuntimeError("please provide the vasp_path for reproduction")
vasp_path = os.path.abspath(self.parameter['vasp_path'])
task_list = reproduce.make_repro(vasp_path, path_to_work)
os.chdir(cwd)
else:
os.chdir(path_to_work)
if os.path.isfile('POSCAR'):
os.remove('POSCAR')
os.symlink(os.path.relpath(equi_contcar), 'POSCAR')
# task_poscar = os.path.join(output, 'POSCAR')
for ii in range(len(all_slabs)):
output_task = os.path.join(path_to_work, 'task.%06d' % ii)
os.makedirs(output_task, exist_ok=True)
os.chdir(output_task)
for jj in ['INCAR', 'POTCAR', 'POSCAR', 'conf.lmp', 'in.lammps']:
if os.path.exists(jj):
os.remove(jj)
task_list.append(output_task)
print("# %03d generate " % ii, output_task, " \t %d atoms" % len(all_slabs[ii].sites))
# make confs
all_slabs[ii].to('POSCAR', 'POSCAR.tmp')
vasp.regulate_poscar('POSCAR.tmp', 'POSCAR')
vasp.sort_poscar('POSCAR', 'POSCAR', ptypes)
vasp.perturb_xz('POSCAR', 'POSCAR', self.pert_xz)
# record miller
np.savetxt('miller.out', all_slabs[ii].miller_index, fmt='%d')
os.chdir(cwd)
return task_list
def task_type(self):
return self.parameter['type']
def task_param(self):
return self.parameter
def _compute_lower(self,
output_file,
all_tasks,
all_res):
output_file = os.path.abspath(output_file)
res_data = {}
ptr_data = os.path.dirname(output_file) + '\n'
if not self.reprod:
ptr_data += "Miller_Indices: \tSurf_E(J/m^2) EpA(eV) equi_EpA(eV)\n"
for ii in all_tasks:
with open(os.path.join(ii, 'inter.json')) as fp:
idata = json.load(fp)
inter_type = idata['type']
equi_path = os.path.abspath(os.path.join(os.path.dirname(output_file), '../relaxation'))
structure_dir = os.path.basename(ii)
if inter_type == 'vasp':
equi_outcar = os.path.join(equi_path, 'OUTCAR')
equi_natoms, equi_epa, equi_vpa = vasp.get_nev(equi_outcar)
outcar = os.path.join(ii, 'OUTCAR')
natoms, epa, vpa = vasp.get_nev(outcar)
if self.static:
e0 = np.array(vasp.get_energies(outcar)) / natoms
epa = e0[0]
boxes = vasp.get_boxes(outcar)
AA = np.linalg.norm(np.cross(boxes[0][0], boxes[0][1]))
elif inter_type in ['deepmd', 'meam', 'eam_fs', 'eam_alloy']:
equi_log = os.path.join(equi_path, 'log.lammps')
equi_natoms, equi_epa, equi_vpa = lammps.get_nev(equi_log)
lmp_log = os.path.join(ii, 'log.lammps')
natoms, epa, vpa = lammps.get_nev(lmp_log)
AA = lammps.get_base_area(lmp_log)
else:
raise RuntimeError('interaction type not supported')
Cf = 1.60217657e-16 / (1e-20 * 2) * 0.001
evac = (epa * natoms - equi_epa * natoms) / AA * Cf
ptr_data += "%s: \t%7.3f %8.3f %8.3f\n" % (structure_dir, evac, epa, equi_epa)
res_data[structure_dir] = [evac, epa, equi_epa]
else:
if 'vasp_path' not in self.parameter:
raise RuntimeError("please provide the vasp_path for reproduction")
vasp_path = os.path.abspath(self.parameter['vasp_path'])
res_data, ptr_data = reproduce.post_repro(vasp_path, all_tasks, ptr_data)
with open(output_file, 'w') as fp:
json.dump(res_data, fp, indent=4)
return res_data, ptr_data