MDStressLab++
Loading...
Searching...
No Matches
readScript.cpp
Go to the documentation of this file.
1/*
2 * main.cpp
3 *
4 * Created on: Nov 3, 2019
5 * Author: Nikhil Admal
6 */
7#include "Mls.h"
8#include "MethodLdad.h"
9#include "MethodSphere.h"
10#include <string>
11#include <iostream>
12#include <tuple>
13#include <cfloat>
14#include <fstream>
15#include "BoxConfiguration.h"
16#include "calculateStress.h"
17#include "Grid.h"
18#include "typedef.h"
19#include <regex>
20
21
22int main()
23{
24 std::string line;
25 auto next_line = [&]() -> std::string {
26 while (std::getline(std::cin, line)) {
27 if (line.empty() || line[0] == '#' || line.find_first_not_of(" \t\r\n") == std::string::npos) continue;
28 return line;
29 }
30 return ""; // EOF
31 };
32
33 // Read LAMMPS data file
34 std::istringstream configFileStream(next_line());
35
36 std::string currentConfigFileName, referenceConfigFileName;
37 std::ifstream referenceConfigFile, currentConfigFile;
38 if (configFileStream >> currentConfigFileName) {
39 std::cout << "LAMMPS data file for current configuration: " << currentConfigFileName << std::endl;
40 currentConfigFile.open(currentConfigFileName);
41 if(!currentConfigFile) MY_ERROR("ERROR: current configuration file could not be opened for reading!");
42 }
43 if (configFileStream >> referenceConfigFileName) {
44 std::cout << "LAMMPS data file for reference configuration: " << referenceConfigFileName << std::endl;
45 referenceConfigFile.open(referenceConfigFileName);
46 if(!referenceConfigFile) MY_ERROR("ERROR: reference configuration file could not be opened for reading!");
47 }
48 if(referenceConfigFileName.empty())
49 referenceConfigFileName= currentConfigFileName;
50
51
52 // get number of particles
53 int numberOfParticles;
54 double xlo, xhi, ylo, yhi, zlo, zhi;
55 while (std::getline(currentConfigFile, line)) {
56 std::string loweredLine = line;
57 std::transform(loweredLine.begin(), loweredLine.end(), loweredLine.begin(), ::tolower);
58
59 if (loweredLine.find("atoms") != std::string::npos && (std::stringstream(line) >> numberOfParticles) ) {
60 }
61 else if (loweredLine.find("xlo xhi") != std::string::npos) {
62 std::istringstream ss(line);
63 ss >> xlo >> xhi;
64 } else if (loweredLine.find("ylo yhi") != std::string::npos) {
65 std::istringstream ss(line);
66 ss >> ylo >> yhi;
67 } else if (loweredLine.find("zlo zhi") != std::string::npos) {
68 std::istringstream ss(line);
69 ss >> zlo >> zhi;
70 break;
71 }
72 }
73
74 if (numberOfParticles < 0) MY_ERROR("Error: Negative number of particles.\n");
75
76 BoxConfiguration body{numberOfParticles, 1};
77 body.readLMP(currentConfigFileName,referenceConfigFileName);
78
79 // Read periodic boundaries
80 Vector3i pbc;
81 std::istringstream pbStream(next_line());
82 pbStream >> pbc(0) >> pbc(1) >> pbc(2);
83 std::cout << "PBCs: " << pbc << std::endl;
84 body.pbc= pbc;
85
86 std::ofstream writeConfigReference("reference.txt");
87 std::ofstream writeConfigCurrent("current.txt");
88 for (int i=0; i<body.numberOfParticles; ++i)
89 {
90 writeConfigReference << body.species[i] << " " << body.coordinates[Reference].row(i) << std::endl;
91 writeConfigCurrent << body.species[i] << " " << body.coordinates[Current].row(i) << std::endl;
92 }
93
94 // Read KIM ID
95 std::string kimID = next_line();
96 std::cout << "KIM ID: " << kimID << std::endl;
97 Kim kim(kimID);
98
99 // Read atom types
100 std::istringstream atomStream(next_line());
101 std::vector<std::string> atomTypes;
102 std::string atom;
103 while (atomStream >> atom) atomTypes.push_back(atom);
104 std::cout << "Atom types: ";
105 for (const auto& a : atomTypes) std::cout << a << " ";
106 std::cout << std::endl;
107
108 // Number of grids
109 int numGrids = std::stoi(next_line());
110 std::cout << "Number of grids: " << numGrids << "\n";
111 Vector3d lowerLimit, upperLimit;
112 int ngridx, ngridy, ngridz;
113 double deltax, deltay, deltaz;
114 std::string outPrefix;
115
116 for (int i = 0; i < numGrids; ++i) {
117 //-------- Read stress method
118 MY_HEADING("Reading stress and grid parameters from input file");
119 bool project= false;
120 std::istringstream stressStream(next_line());
121 std::string stressMethod, averagingDomain, averagingDomainParameter;
122 double averagingDomainSize;
123
124 stressStream >> stressMethod >> averagingDomain ;
125 std::cout << "stress method: " << stressMethod << "\n";
126 std::cout << "averaging domain: " << averagingDomain << "\n";
127 if (stressMethod == "project")
128 project= true;
129
130 Vector3d xdir, ydir, zdir;
131 Matrix3d ldadVectors;
132
133 if (averagingDomain=="ldad")
134 {
135 stressStream >> averagingDomainParameter;
136 std::cout << "averaging domain parameter: " << averagingDomainParameter << "\n";
137
138 // read crystallographic directions and domain size to calculate ldadVectors
139 if (averagingDomainParameter=="bcc" || averagingDomainParameter=="fcc")
140 {
141 stressStream >> averagingDomainSize >> xdir(0) >> xdir(1) >> xdir(2) >>
142 ydir(0) >> ydir(1) >> ydir(2) >>
143 zdir(0) >> zdir(1) >> zdir(2);
144 std::cout << "averaging domain size : " << averagingDomainSize << "\n";
145
146 Vector3d basis1, basis2, basis3;
147 if (averagingDomainParameter=="bcc"){
148 basis1 << -0.5,0.5,0.5; basis2 << 0.5,-0.5,0.5; basis3 << 0.5,0.5,-0.5;
149 }
150 else if (averagingDomainParameter== "fcc"){
151 basis1 << 0,0.5,0.5; basis2 << 0.5,0,0.5; basis3 << 0.5,0.5,0;
152 }
153 Matrix3d rotation;
154
155 // calculate lattice vectors
156 rotation.row(0)= xdir.normalized();
157 rotation.row(1)= ydir.normalized();
158 rotation.row(2)= zdir.normalized();
159 assert((rotation.transpose()*rotation).isIdentity());
160 ldadVectors.col(0)= rotation*basis1.transpose();
161 ldadVectors.col(1)= rotation*basis2.transpose();
162 ldadVectors.col(2)= rotation*basis3.transpose();
163 ldadVectors= ldadVectors* averagingDomainSize;
164 }
165 // instead directly read ldad lattice vectors
166 else if (averagingDomainParameter=="lat")
167 stressStream >> ldadVectors(0,0) >> ldadVectors(1,0) >> ldadVectors(2,0) >>
168 ldadVectors(0,1) >> ldadVectors(1,1) >> ldadVectors(2,1) >>
169 ldadVectors(0,2) >> ldadVectors(1,2) >> ldadVectors(2,2);
170 }
171 else if (averagingDomain=="sphere"){
172 stressStream >> averagingDomainSize;
173 std::cout << "averaging domain size : " << averagingDomainSize << "\n";
174 }
175
176 //------------ Read grid
177 std::istringstream gridStream(next_line());
178 std::string token[9];
179 for (auto & j : token) {
180 gridStream >> j;
181 }
182 // Lambda to convert token or use fallback
183 auto parseOrFallback = [](const std::string& s, double fallback) -> double {
184 return (s == "*") ? fallback : std::stod(s);
185 };
186 // Use fallback values if needed
187 lowerLimit[0]= parseOrFallback(token[0], xlo)-FLT_EPSILON;
188 lowerLimit[1]= parseOrFallback(token[1], ylo)-FLT_EPSILON;
189 lowerLimit[2]= parseOrFallback(token[2], zlo)-FLT_EPSILON;
190 //lowerLimit[0]= parseOrFallback(token[0], xlo);
191 //lowerLimit[1]= parseOrFallback(token[1], ylo);
192 //lowerLimit[2]= parseOrFallback(token[2], zlo);
193
194 upperLimit[0]= parseOrFallback(token[3], xhi);
195 upperLimit[1]= parseOrFallback(token[4], yhi);
196 upperLimit[2]= parseOrFallback(token[5], zhi);
197
198 // Parse delta values and diameter directly
199 deltax = std::stod(token[6]);
200 deltay = std::stod(token[7]);
201 deltaz = std::stod(token[8]);
202
203 ngridx= (abs(deltax) > FLT_EPSILON) ? floor((upperLimit(0)-lowerLimit(0))/deltax) : 1;
204 ngridy= (abs(deltay) > FLT_EPSILON) ? floor((upperLimit(1)-lowerLimit(1))/deltay) : 1;
205 ngridz= (abs(deltaz) > FLT_EPSILON) ? floor((upperLimit(2)-lowerLimit(2))/deltaz) : 1;
206 //ngridx= (abs(deltax) > FLT_EPSILON) ? floor((upperLimit(0)-lowerLimit(0)+FLT_EPSILON)/deltax) : 1;
207 //ngridy= (abs(deltay) > FLT_EPSILON) ? floor((upperLimit(1)-lowerLimit(1)+FLT_EPSILON)/deltay) : 1;
208 //ngridz= (abs(deltaz) > FLT_EPSILON) ? floor((upperLimit(2)-lowerLimit(2)+FLT_EPSILON)/deltaz) : 1;
209
210 std::cout << "Grid Limits: ";
211 std::cout << lowerLimit << " " << upperLimit << std::endl;
212 std::cout << "Number of grid points: " << ngridx << " " << ngridy << " " << ngridz << std::endl;
213
214 // ----------- Output prefix
215 outPrefix= next_line();
216 std::cout << "Output file: " << outPrefix << "\n";
217
218
219 if (averagingDomain=="ldad") {
220 try {
221 Grid<Reference> grid(lowerLimit, upperLimit, ngridx, ngridy, ngridz);
222 //Matrix3d ldadVectors= averagingDomainSize*Matrix3d::Identity();
223
224 MethodLdadTrigonometric ldadDomain(ldadVectors);
225 Stress<MethodLdadTrigonometric, Piola> ldadTrigonometricStress(ldadDomain, &grid);
226 calculateStress(body, kim,
227 std::tie(ldadTrigonometricStress),
228 std::tie(), project);
229
230 double mlsRadius= 10.0;
231 Mls mls(body,&grid,mlsRadius,outPrefix);
232 std::vector<Matrix3d> cauchyPushedField;
233 mls.pushToCauchy(ldadTrigonometricStress.field,cauchyPushedField);
234 mls.writePushedCauchyStress(cauchyPushedField);
236 }
237 catch (const std::runtime_error &e) {
238 std::cout << e.what() << std::endl;
239 std::cout << "Compute stress with projected forces failed. Moving on" << std::endl;
240 }
241 }
242 else if(averagingDomain=="sphere"){
243 Grid<Current> grid(lowerLimit, upperLimit, ngridx, ngridy, ngridz);
244 MethodSphere hardy(averagingDomainSize, "hardy");
245
246 // stress calculation using projected forces
247 try {
248 Stress<MethodSphere, Cauchy> hardyStress(hardy, &grid);
249
250 calculateStress(body, kim,
251 std::tie(),
252 std::tie(hardyStress), project);
253 hardyStress.write(outPrefix);
254 }
255 catch (const std::runtime_error &e) {
256 std::cout << e.what() << std::endl;
257 std::cout << "Compute stress with projected forces failed. Moving on" << std::endl;
258 }
259 }
260 else
261 MY_ERROR("Unknown stress method: " + stressMethod);
262 }
263
264 return 0;
265}
int calculateStress(const BoxConfiguration &body, Kim &kim, std::tuple<> stress, const bool &projectForces=false)
Represents a particle configuration including simulation box information.
void readLMP(const std::string &, const ConfigType &configType)
Reads a configuration from an LAMMPS-style dump file.
Definition Grid.h:31
Definition kim.h:21
Implements radially symmetric weighting functions (Hardy, Virial) and its associated bond function fo...
Computes deformation gradients and stress tensors using Moving Least Squares (MLS) interpolation on a...
Definition Mls.h:26
void writeDeformationGradient()
Writes the deformation gradient at each grid point to a file.
Definition Mls.cpp:832
void pushToCauchy(const std::vector< Matrix3d > &piolaStress, std::vector< Matrix3d > &cauchyStress)
Converts Piola-Kirchhoff stress to Cauchy stress at grid points.
Definition Mls.cpp:796
void writePushedCauchyStress(std::vector< Matrix3d > &cauchyStress)
Writes Cauchy stresses at grid points to a file.
Definition Mls.cpp:887
std::vector< Matrix3d > field
A three-dimensional stress field.
Definition Stress.h:39
void write()
Definition Stress.h:90
int main()
#define MY_ERROR(message)
Definition typedef.h:17
Eigen::Matrix< double, DIM, DIM, Eigen::RowMajor > Matrix3d
Definition typedef.h:56
Eigen::Matrix< double, 1, DIM, Eigen::RowMajor > Vector3d
Definition typedef.h:60
Eigen::Matrix< int, 1, DIM, Eigen::RowMajor > Vector3i
Definition typedef.h:61
@ Current
Definition typedef.h:70
@ Reference
Definition typedef.h:69
#define MY_HEADING(heading)
Definition typedef.h:36