oILAB
Loading...
Searching...
No Matches
Lammps.h
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 11/19/24.
3//
4
5#ifndef OILAB_LAMMPS_H
6#define OILAB_LAMMPS_H
7
8#include <iostream>
9#include <fstream>
10#include <vector>
11#include <string>
12#include <cstdlib>
13#include <cmath>
14#include <sstream>
15#include <sys/stat.h>
16#include <unistd.h>
17#include <Eigen/Eigen>
18#include <iomanip>
19
20std::tuple<Eigen::MatrixXd,
21 Eigen::Matrix3d,
22 Eigen::Vector3d> read_oILAB_output(const std::string& path)
23{
24 Eigen::MatrixXd atoms;
25 Eigen::Matrix3d box;
26 box.setZero();
27 Eigen::Vector<double,9> boxVectorized;
28 boxVectorized.setZero();
29 Eigen::Vector3d origin;
30 origin.setZero();
31
32 std::ifstream file(path);
33 if (!file.is_open()) {
34 std::cerr << "Error opening file for reading oilab config file: " << path << std::endl;
35 }
36
37 std::string line;
38 int number_atoms = 0;
39
40 int lineCount= 0;
41 int atomCount= 0;
42 int boxCount= -1;
43 int originCount= -1;
44
45 while (std::getline(file, line)) {
46 lineCount++;
47 std::stringstream ss(line);
48 std::string field;
49 std::vector<std::string> fields;
50
51 while (ss >> field) {
52 fields.push_back(field);
53 }
54
55 if (lineCount == 1) {
56 number_atoms = std::stoi(fields[0]);
57 atoms.resize(number_atoms,5);
58 atoms.setZero();
59 }
60 else if (lineCount == 2)
61 {
62 for (size_t i = 0; i < fields.size(); ++i) {
63 if (fields[i].find("Lattice=") != std::string::npos) {
64 boxCount++;
65 continue;
66 }
67 if (boxCount>=0 && boxCount<9) {
68 boxVectorized(boxCount) = std::stod(fields[i]);
69 boxCount++;
70 continue;
71 }
72 if (fields[i].find("origin=") != std::string::npos)
73 {
74 originCount++;
75 continue;
76 }
77 if (originCount>-1 && originCount<3) {
78 origin(originCount) = std::stod(fields[i]);
79 originCount++;
80 continue;
81 }
82 }
83 }
84 else if (lineCount > 2) {
85 // Reading atom data (type, x, y, z, radius)
86 Eigen::VectorXd atom_data(5);
87 int coordCount= 0;
88 for (const auto &field : fields) {
89 atom_data(coordCount)= std::stod(field);
90 coordCount++;
91 }
92 atoms.row(atomCount)= atom_data;
93 atomCount++;
94 }
95 }
96 box= boxVectorized.reshaped(3,3);
97
98 return {atoms, box, origin};
99
100}
101
102/*-----------------------*/
103//std::string write_lammps_datafile(const std::string &filename, const std::vector<std::vector<double>> &box, const std::vector<std::vector<double>> &atom_data, int atom_types) {
104void write_lammps_datafile(const std::string &filename,
105 const std::vector<std::vector<double>> &box,
106 const Eigen::MatrixXd& atom_data,
107 int atom_types)
108{
109 std::ofstream file(filename);
110 if (!file.is_open())
111 std::cerr << "Error opening file for writing lammps configuration file: " << filename << std::endl;
112
113 file << "# LAMMPS data file via write_data\n";
114 file << "\n";
115 file << atom_data.rows() << " atoms\n";
116 file << atom_types << " atom types\n";
117 file << "\n";
118 file << std::setprecision(8) << box[0][0] << " " << std::setprecision(8) << box[0][1] << " xlo xhi\n";
119 file << std::setprecision(8) << box[1][0] << " " << std::setprecision(8) << box[1][1] << " ylo yhi\n";
120 file << std::setprecision(8) << box[2][0] << " " << std::setprecision(8) << box[2][1] << " zlo zhi\n";
121 file << "0 0 0 xy xz yz\n";
122 file << "\n";
123 file << "Atoms # atomic\n";
124 file << "\n";
125 for (size_t i = 0; i < atom_data.rows(); ++i) {
126 //file << atom_data(i,0) << " " << atom_data(i,1) << " " << atom_data(i,2) << " " << atom_data(i,3) << " " << atom_data(i,4) << " 0 0 0\n";
127 file << std::setprecision(8) << atom_data.row(i) << " 0 0 0\n";
128 }
129}
130
131/*-----------------------*/
132void write_lammps_input_script(const std::string &filename,
133 const std::string &infile,
134 const std::string &outfile,
135 double gb_thickness_parameter,
136 const std::string &potential_file_path,
137 const std::string &output_dump_file) {
138 std::ofstream file(filename);
139 if (!file.is_open()) {
140 std::cerr << "Error opening file for writing lammps input script: " << filename << std::endl;
141 return;
142 }
143
144 file << "# Find energy of a config and write in " << outfile << " file\n";
145 file << "\n";
146 file << "variable potential_path string " << potential_file_path << "\n";
147 file << "clear\n";
148 file << "units metal\n";
149 file << "boundary p p p\n";
150 file << "atom_style atomic\n";
151 file << "neighbor 1.0 bin\n";
152 file << "neigh_modify every 1 delay 2 check yes\n";
153 file << "read_data " << infile << "\n";
154 file << "pair_style eam/alloy\n";
155 file << "pair_coeff * * ${potential_path} Cu Cu\n";
156 file << "delete_atoms overlap 1e-2 all all\n";
157 file << "variable area equal ly*lz\n";
158 file << "variable xlogb equal xlo+" << std::setprecision(8) << gb_thickness_parameter << "\n";
159 file << "variable xhigb equal xhi-" << std::setprecision(8) << gb_thickness_parameter << "\n";
160 file << "variable xlobulk equal xlo+" << std::setprecision(8) << gb_thickness_parameter << "\n";
161 file << "variable xhibulk equal xlo+" << std::setprecision(8) << 1.5*gb_thickness_parameter << "\n";
162 file << "region GB block ${xlogb} ${xhigb} INF INF INF INF side in units box\n";
163 file << "region BULK block ${xlobulk} ${xhibulk} INF INF INF INF side in units box\n";
164 file << "group GB region GB\n";
165 file << "group BULK region BULK\n";
166 file << "compute peratom GB pe/atom\n";
167 file << "compute peratombulk BULK pe/atom\n";
168 file << "compute pe GB reduce sum c_peratom\n";
169 file << "compute pebulk GB reduce sum c_peratombulk\n";
170 file << "variable peGB equal c_pe\n";
171 file << "variable peBULK equal c_pebulk\n";
172 file << "variable atomsGB equal count(GB)\n";
173 file << "variable atomsBULK equal count(BULK)\n";
174 file << "thermo 1000\n";
175 file << "compute 1 all ke/atom\n";
176 file << "compute cna all cna/atom 3.08133\n";
177 file << "compute csys all centro/atom fcc\n";
178 file << "compute 3 all pe/atom\n";
179 file << "compute 4 all stress/atom NULL pair\n";
180 file << "timestep 0.001\n";
181 file << "thermo_style custom step temp ke pe etotal press pxx pyy pzz pxy pxz pyz ly lx lz xy xz yz c_pe v_atomsGB v_peBULK v_atomsBULK\n";
182 file << "dump OUT0 all custom 10 " << output_dump_file << " id type x y z fx fy fz c_3 c_1 vx vy vz c_4[1] c_4[2] c_4[3] c_4[4] c_4[5] c_4[6]\n";
183 file << "run 0\n";
184 file << "variable coh equal (${peBULK}/${atomsBULK})\n";
185 file << "variable GBene equal (${peGB}-${coh}*${atomsGB})\n";
186 file << "print \"coh = ${coh} energy = ${peGB} numAtoms = ${atomsGB} GBene = ${GBene} area = ${area}\" file " << outfile << "\n";
187 file << "\n";
188}
189
190/*-----------------------*/
191std::vector<std::vector<double>> read_python_outfile(const std::string &path) {
192 std::vector<std::vector<double>> data;
193 std::ifstream file(path);
194 if (!file.is_open()) {
195 std::cerr << "Error opening file for reading lammps output: " << path << std::endl;
196 return data;
197 }
198
199 std::string line;
200 while (std::getline(file, line)) {
201 std::stringstream ss(line);
202 std::string field;
203 double state_id, area, total_energy, gb_energy, gb_density;
204
205 while (ss >> field)
206 {
207 if (field == "coh") {
208 ss >> field;
209 ss >> state_id;
210 }
211 else if (field == "energy") {
212 ss >> field;
213 ss >> total_energy;
214 }
215 else if (field == "GBene") {
216 ss >> field;
217 ss >> gb_energy;
218 }
219 else if (field == "numAtoms") {
220 ss >> field;
221 ss >> gb_density;
222 }
223 else if (field == "area")
224 {
225 ss >> field;
226 ss >> area;
227 break;
228 }
229 }
230
231 //gb_density = gb_density / area;
232 data.push_back({state_id, area, gb_energy, gb_density});
233 }
234
235 return data;
236}
237
238
239std::pair<double, double> energy(const std::string& lammpsLocation,
240 const std::string& oilabConfigFile,
241 const std::string& potentialFile)
242{
243 // Write data
244 std::string threadNumber= std::to_string(omp_get_thread_num());
245 std::string lammpsInputFile= "in"+ threadNumber +".find_energy";
246 std::string lammpsDataFile= "data" + threadNumber + ".lammps_input";
247 std::string lammpsDumpFile= "dump" + threadNumber + ".lammpsConfigs";
248 std::string outfile = "lmp_mesostate_energies" + threadNumber + ".txt";
249
250 // Read data
251 auto [atoms, box, origin] = read_oILAB_output(oilabConfigFile);
252
253 // Find rotation and new box
254 Eigen::Matrix3d R= box.transpose();
255 for(int i=0; i<3; ++i)
256 R.row(i).normalize();
257 Eigen::Matrix3d new_box= R*box;
258
259 // New atoms
260 Eigen::MatrixXd new_atoms(atoms.rows(),5);
261 for (size_t i = 0; i < atoms.rows(); ++i) {
262 new_atoms(i,Eigen::seq(2,4))= R*(atoms(i,Eigen::seq(1,3)).transpose());
263 /*
264 new_atoms[i][2] = new_atoms[i][2] + R[0][0]*atoms[i][1] + R[0][1]*atoms[i][2] + R[0][2]*atoms[i][3];
265 new_atoms[i][3] = new_atoms[i][3] + R[1][0]*atoms[i][1] + R[1][1]*atoms[i][2] + R[1][2]*atoms[i][3];
266 new_atoms[i][4] = new_atoms[i][4] + R[2][0]*atoms[i][1] + R[2][1]*atoms[i][2] + R[2][2]*atoms[i][3];
267 */
268 new_atoms(i,0)= i + 1;
269 new_atoms(i,1) = atoms(i,0);
270 }
271
272 // Determine nbox values
273 std::vector<std::vector<double>> nbox(3, std::vector<double>(2, 0.0));
274 for (size_t i = 0; i < 3; ++i) {
275 if (i == 0) {
276 nbox[i][0] = -new_box(i,i) / 2;
277 nbox[i][1] = new_box(i,i) / 2;
278 } else if (i == 1) {
279 //nbox[i][0] = origin(i);
280 nbox[i][0] = 0.0;
281 nbox[i][1] = new_box(i,i);
282 } else {
283 nbox[i][0] = 0.0;
284 nbox[i][1] = new_box(i,i);
285 }
286 }
287
288 //double gb_thickness_parameter = 0.75 * (new_box(0,0) * 0.5);
289 double gb_thickness_parameter= 6;
290
291 // Write files
292 write_lammps_datafile(lammpsDataFile, nbox, new_atoms, 2);
293 write_lammps_input_script(lammpsInputFile, lammpsDataFile, outfile, gb_thickness_parameter, potentialFile, lammpsDumpFile);
294
295 // Run the LAMMPS script
296 std::string command = lammpsLocation +" -in " + lammpsInputFile + " > /dev/null 2>&1";
297 std::system(command.c_str());
298
299 // Read energy
300 auto data_energy = read_python_outfile(outfile);
301
302 return {data_energy[0][3], data_energy[0][2]};
303}
304
305#endif //OILAB_LAMMPS_H
std::vector< std::vector< double > > read_python_outfile(const std::string &path)
Definition Lammps.h:191
std::tuple< Eigen::MatrixXd, Eigen::Matrix3d, Eigen::Vector3d > read_oILAB_output(const std::string &path)
Definition Lammps.h:22
std::pair< double, double > energy(const std::string &lammpsLocation, const std::string &oilabConfigFile, const std::string &potentialFile)
Definition Lammps.h:239
void write_lammps_input_script(const std::string &filename, const std::string &infile, const std::string &outfile, double gb_thickness_parameter, const std::string &potential_file_path, const std::string &output_dump_file)
Definition Lammps.h:132
void write_lammps_datafile(const std::string &filename, const std::vector< std::vector< double > > &box, const Eigen::MatrixXd &atom_data, int atom_types)
Definition Lammps.h:104