-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathphonopy.cpp
334 lines (294 loc) · 10.3 KB
/
phonopy.cpp
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
#ifdef FFTW3
#include "phonopy.h"
#include "global.h"
#include "dynmat.h"
#include "input.h"
#include "memory.h"
#include <fftw3.h>
#include <cmath>
#include <cstdio>
#include <cstdlib>
#include <cstring>
#include <map>
/* ----------------------------------------------------------------------------
* Class Phonopy is designed to interface with phonopy.
* ---------------------------------------------------------------------------- */
Phonopy::Phonopy(DynMat *dynmat)
{
dm = dynmat;
memory = new Memory();
sysdim = dm->sysdim;
fftdim = dm->fftdim;
input = dm->input;
fftdim2 = fftdim * fftdim;
nucell = dm->nucell;
nx = ny = nz = 5;
write(1);
char str[MAXLINE];
if (input == NULL) input = new UserInput(0);
input->read_stdin(str);
if (count_words(str) >= 3){
nx = atoi(strtok(str," \t\n\r\f"));
ny = atoi(strtok(NULL," \t\n\r\f"));
nz = atoi(strtok(NULL," \t\n\r\f"));
}
if (dm->nx == 1) nx = 1;
if (dm->ny == 1) ny = 1;
if (dm->nz == 1) nz = 1;
if (nx < 1) nx = 1;
if (ny < 1) ny = 1;
if (nz < 1) nz = 1;
npt = nx * ny * nz;
write(2);
memory->create(mass, nucell, "Phonopy:mass");
for (int i = 0; i < nucell; ++i){
double m = 1.0/dm->M_inv_sqrt[i];
mass[i] = m * m;
}
memory->create(FC_all, npt, fftdim2, "phonopy:FC_all");
get_my_FC();
phonopy();
return;
}
/* ----------------------------------------------------------------------------
* Deconstructor, free the memory used.
* ---------------------------------------------------------------------------- */
Phonopy::~Phonopy()
{
memory->destroy(mass);
memory->destroy(FC_all);
delete memory;
dm = NULL;
return;
}
/* ----------------------------------------------------------------------------
* Subroutine to write the outputs to screen.
* ---------------------------------------------------------------------------- */
void Phonopy::write(int flag)
{
if (flag == 1){ // basic information
puts("================================================================================");
printf("Now to prepare the input files for phonopy.\n");
printf("The dimension of your present supercell is : %d x %d x %d.\n", dm->nx, dm->ny, dm->nz);
printf("The size of the force constant matrix will be: %d x %d.\n", dm->npt*3, dm->npt*3);
printf("Please define your desired dimension [5 5 5] : ");
} else if (flag == 2){
printf("\nThe new dimension of the supercell will be : %d x %d x %d.\n", nx, ny, nz);
printf("The size of the force constant matrix is then: %d x %d.\n", npt*3, npt*3);
} else if (flag == 3){
printf("\nNow to prepare the phonopy FORCE_CONSTANTS ..."); fflush(stdout);
} else if (flag == 4){
printf("Done!\nThe force constants information is extracted and written to FORCE_CONSTANTS,\n");
printf("the primitive cell is written to POSCAR.primitive, and the input file for\n");
printf("phonopy band evaluation is written to band.conf.\n\n");
printf("One should be able to obtain the phonon band structure after\n");
printf(" 1) Correcting the `element names` in POSCAR.primitive and band.conf;\n");
printf(" 2) Running `phonopy --readfc -c POSCAR.primitive -p band.conf`.\n\n");
printf("Or the phonon density of states after\n");
printf(" 1) Correcting the `element names` in POSCAR.primitive and mesh.conf;\n");
printf(" 2) Running `phonopy --readfc -c POSCAR.primitive -p mesh.conf`.\n");
puts("--------------------------------------------------------------------------------");
printf("*** Remember to modify the `element names`. ***\n");
} else if (flag == 5){
puts("================================================================================");
}
return;
}
/* ----------------------------------------------------------------------------
* Driver to obtain the force constant matrix
* ---------------------------------------------------------------------------- */
void Phonopy::get_my_FC()
{
double q[3];
int ipt = 0;
for (int ix = 0; ix < nx; ++ix)
for (int iy = 0; iy < ny; ++iy)
for (int iz = 0; iz < nz; ++iz){
q[0] = double(ix)/double(nx);
q[1] = double(iy)/double(ny);
q[2] = double(iz)/double(nz);
dm->getDMq(q);
int ndim = 0;
for (int idim = 0; idim < fftdim; ++idim)
for (int jdim = 0; jdim < fftdim; ++jdim){
FC_all[ipt][ndim] = dm->DM_q[idim][jdim];
double m = sqrt(mass[idim/sysdim] * mass[jdim/sysdim]);
FC_all[ipt][ndim].r *= m;
FC_all[ipt][ndim].i *= m;
++ndim;
}
++ipt;
}
return;
}
/* ----------------------------------------------------------------------------
* Method to write out the force constants and related files for
* postprocessing with phonopy.
* ---------------------------------------------------------------------------- */
void Phonopy::phonopy()
{
// output info
write(3);
fftw_complex *in, *out;
double **fc;
memory->create(in, npt, "phonopy:in");
memory->create(out, npt, "phonopy:in");
memory->create(fc, npt, fftdim2, "phonopy:in");
fftw_plan plan = fftw_plan_dft_3d(nx, ny, nz, in, out, -1, FFTW_ESTIMATE);
double factor = dm->eml2fc / double(npt);
for (int idim = 0; idim < fftdim2; ++idim){
for (int i = 0; i < npt; ++i){
in[i][0] = FC_all[i][idim].r;
in[i][1] = FC_all[i][idim].i;
}
fftw_execute(plan);
for (int i = 0; i < npt; ++i) fc[i][idim] = out[i][0] * factor;
}
fftw_destroy_plan(plan);
memory->destroy(in);
memory->destroy(out);
// in POSCAR, atoms are sorted/aggregated by type, while for LAMMPS there is no such requirment
int *type_id = new int[nucell];
int *num_type = new int[nucell];
int ntype = 0;
for (int i = 0; i < nucell; ++i) num_type[i] = 0;
for (int i = 0; i < nucell; ++i){
int ip = ntype;
for (int j = 0; j < ntype; ++j){
if (dm->attyp[i] == type_id[j]) ip = j;
}
if (ip == ntype) type_id[ntype++] = dm->attyp[i];
num_type[ip]++;
}
std::map<int, int> iu_by_type;
iu_by_type.clear();
int id_new = 0;
for (int i = 0; i < ntype; ++i){
for (int j = 0; j < nucell; ++j){
if (dm->attyp[j] == type_id[i]) iu_by_type[id_new++] = j;
}
}
// write the FORCE_CONSTANTS file
FILE *fp = fopen("FORCE_CONSTANTS", "w");
int natom = npt * nucell;
fprintf(fp, "%d %d\n", natom, natom);
for (int i = 0; i < natom; ++i){
int iu = i / npt;
int iz = (i % npt) / (nx * ny);
int iy = (i % (nx *ny) ) / nx;
int ix = i % nx;
iu = iu_by_type[iu];
for (int j = 0; j < natom; ++j){
int ju = j / npt;
int jz = (j % npt) / (nx * ny);
int jy = (j % (nx *ny) ) / nx;
int jx = j % nx;
int dx = abs(ix - jx);
int dy = abs(iy - jy);
int dz = abs(iz - jz);
int id = (dx * ny + dy) * nz + dz;
ju = iu_by_type[ju];
fprintf(fp, "%d %d\n", i+1, j+1);
for (int idim = iu * sysdim; idim < (iu+1)*sysdim; ++idim){
for (int jdim = ju * sysdim; jdim < (ju+1)*sysdim; ++jdim){
int dd = idim * fftdim + jdim;
fprintf(fp, " %lg", fc[id][dd]);
}
fprintf(fp, "\n");
}
}
}
fclose(fp);
iu_by_type.clear();
memory->destroy(fc);
// write the primitive cell in POSCAR format
fp = fopen("POSCAR.primitive", "w");
fprintf(fp, "Fix-phonon unit cell");
for (int ip = 0; ip < ntype; ++ip){
for (int i = 0; i < nucell; ++i){
if (dm->attyp[i] == type_id[ip]){
fprintf(fp, ", Elem-%d: %lg", type_id[ip], mass[i]);
break;
}
}
}
fprintf(fp, "\n1.\n");
int ndim = 0;
for (int idim = 0; idim < 3; ++idim){
for (int jdim = 0; jdim < 3; ++jdim) fprintf(fp, "%lg ", dm->basevec[ndim++]);
fprintf(fp, "\n");
}
for (int ip = 0; ip < ntype; ++ip) fprintf(fp, "Elem-%d ", type_id[ip]);
fprintf(fp, "\n");
for (int ip = 0; ip < ntype; ++ip) fprintf(fp, "%d ", num_type[ip]);
fprintf(fp, "\nDirect\n");
for (int ip = 0; ip < ntype; ++ip){
for (int i = 0; i < nucell; ++i){
if (dm->attyp[i] == type_id[ip]){
fprintf(fp, "%lg %lg %lg\n", dm->basis[i][0], dm->basis[i][1], dm->basis[i][2]);
}
}
}
fclose(fp);
// mesh.conf
fp = fopen("mesh.conf", "w");
fprintf(fp, "# From Fix-phonon");
for (int ip = 0; ip < ntype; ++ip){
for (int i = 0; i < nucell; ++i){
if (dm->attyp[i] == type_id[ip]){
fprintf(fp, ", Elem-%d: %lg", type_id[ip], mass[i]);
break;
}
}
}
fprintf(fp, "\n\nATOM_NAME = ");
for (int ip = 0; ip < ntype; ++ip) fprintf(fp, "Elem-%d ", type_id[ip]);
fprintf(fp, "\nDIM = %d %d %d\n", nx, ny, nz);
fprintf(fp, "MP = 31 31 31\nFORCE_CONSTANTS = READ\n");
fprintf(fp, "#FC_SYMMETRY = .TRUE.\n#SYMMETRY_TOLERANCE = 0.01\n");
fclose(fp);
// band.conf
fp = fopen("band.conf", "w");
fprintf(fp, "# From Fix-phonon");
for (int ip = 0; ip < ntype; ++ip){
for (int i = 0; i < nucell; ++i){
if (dm->attyp[i] == type_id[ip]){
fprintf(fp, ", Elem-%d: %lg", type_id[ip], mass[i]);
break;
}
}
}
fprintf(fp, "\n\nATOM_NAME = ");
for (int ip = 0; ip < ntype; ++ip) fprintf(fp, "Elem-%d ", type_id[ip]);
fprintf(fp, "\nDIM = %d %d %d\nBAND = AUTO\n", nx, ny, nz);
fprintf(fp, "BAND_POINTS = 21\nFORCE_CONSTANTS = READ\nBAND_CONNECTION = .TRUE.\n");
fprintf(fp, "#FC_SYMMETRY = .TRUE.\n#SYMMETRY_TOLERANCE = 0.01\n");
// output info
write(4);
write(5);
delete[] type_id;
delete[] num_type;
return;
}
/*------------------------------------------------------------------------------
* Method to count # of words in a string, without destroying the string
*----------------------------------------------------------------------------*/
int Phonopy::count_words(const char *line)
{
int n = strlen(line) + 1;
char *copy;
memory->create(copy, n, "count_words:copy");
strcpy(copy,line);
char *ptr;
if ((ptr = strchr(copy,'#'))) *ptr = '\0';
if (strtok(copy," \t\n\r\f") == NULL) {
memory->destroy(copy);
return 0;
}
n = 1;
while (strtok(NULL," \t\n\r\f")) n++;
memory->destroy(copy);
return n;
}
/*----------------------------------------------------------------------------*/
#endif