forked from etmc/tmLQCD
-
Notifications
You must be signed in to change notification settings - Fork 1
/
gauge_io.c
396 lines (372 loc) · 12.6 KB
/
gauge_io.c
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
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
/***********************************************************************
* Copyright (C) 2002,2003,2004,2005,2006,2007,2008 Carsten Urbach
*
* This file is part of tmLQCD.
*
* tmLQCD is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* tmLQCD is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with tmLQCD. If not, see <http://www.gnu.org/licenses/>.
***********************************************************************/
#define _FILE_OFFSET_BITS 64
#include"lime.h"
#ifdef HAVE_CONFIG_H
# include<config.h>
#endif
#include<stdlib.h>
#include<stdio.h>
#include<string.h>
#include<time.h>
#include<sys/time.h>
#include<sys/types.h>
#ifdef MPI
# include<mpi.h>
# include<unistd.h>
#endif
#include<math.h>
#include"global.h"
#include"su3.h"
#include"lime.h"
#include"read_input.h"
#include"io_utils.h"
#include"dml.h"
#include"io.h"
#include"gauge_io.h"
/* #define MAXBUF 1048576 */
int write_binary_gauge_data(LimeWriter * limewriter,
const int prec, DML_Checksum * ans) {
int x, X, y, Y, z, Z, tt, t0, tag=0, id=0, status=0;
su3 tmp[4];
su3 tmp3[4];
float tmp2[72];
int coords[4];
/* n_uint64_t bytes; */
n_uint64_t bytes;
DML_SiteRank rank;
#ifdef MPI
MPI_Status mpi_status;
#endif
DML_checksum_init(ans);
if(prec == 32) bytes = (n_uint64_t)2*sizeof(su3);
else bytes = (n_uint64_t)4*sizeof(su3);
for(t0 = 0; t0 < T*g_nproc_t; t0++) {
tt = t0 - g_proc_coords[0]*T;
coords[0] = t0 / T;
for(z = 0; z < LZ*g_nproc_z; z++) {
Z = z - g_proc_coords[3]*LZ;
coords[3] = z / LZ;
for(y = 0; y < LY*g_nproc_y; y++) {
tag = 0;
Y = y - g_proc_coords[2]*LY;
coords[2] = y / LY;
for(x = 0; x < LX*g_nproc_x; x++) {
X = x - g_proc_coords[1]*LX;
coords[1] = x / LX;
#ifdef MPI
MPI_Cart_rank(g_cart_grid, coords, &id);
#endif
if(g_cart_id == 0) {
/* Rank should be computed by proc 0 only */
rank = (DML_SiteRank) (((t0*LZ*g_nproc_z + z)*LY*g_nproc_y + y)*LX*g_nproc_x + x);
if(g_cart_id == id) {
memcpy(&tmp3[0], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][1], sizeof(su3));
memcpy(&tmp3[1], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][2], sizeof(su3));
memcpy(&tmp3[2], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][3], sizeof(su3));
memcpy(&tmp3[3], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][0], sizeof(su3));
#ifndef WORDS_BIGENDIAN
if(prec == 32) {
byte_swap_assign_double2single(tmp2, tmp3, 4*sizeof(su3)/8);
DML_checksum_accum(ans, rank, (char*) tmp2, 4*sizeof(su3)/2);
status = limeWriteRecordData((void*)&tmp2, &bytes, limewriter);
}
else {
byte_swap_assign(tmp, tmp3, 4*sizeof(su3)/8);
DML_checksum_accum(ans, rank, (char*) tmp, 4*sizeof(su3));
status = limeWriteRecordData((void*)&tmp, &bytes, limewriter);
}
#else
if(prec == 32) {
double2single(tmp2, tmp3, 4*sizeof(su3)/8);
DML_checksum_accum(ans, rank, (char*) tmp2, 4*sizeof(su3)/2);
status = limeWriteRecordData((void*)&tmp2, &bytes, limewriter);
}
else {
DML_checksum_accum(ans, rank, (char*) tmp3, 4*sizeof(su3));
status = limeWriteRecordData((void*)&tmp3, &bytes, limewriter);
}
#endif
}
#ifdef MPI
else {
if(prec == 32) {
MPI_Recv(tmp2, 4*sizeof(su3)/8, MPI_FLOAT, id, tag, g_cart_grid, &mpi_status);
DML_checksum_accum(ans, rank, (char*) tmp2, 4*sizeof(su3)/2);
status = limeWriteRecordData((void*)&tmp2, &bytes, limewriter);
}
else {
MPI_Recv(tmp, 4*sizeof(su3)/8, MPI_DOUBLE, id, tag, g_cart_grid, &mpi_status);
DML_checksum_accum(ans, rank, (char*) tmp, 4*sizeof(su3));
status = limeWriteRecordData((void*)&tmp, &bytes, limewriter);
}
}
#endif
if(status < 0 ) {
fprintf(stderr, "LIME write error %d\n", status);
fprintf(stderr, "x %d, y %d, z %d, t %d (%d,%d,%d,%d)\n",x,y,z,tt,X,Y,Z,tt);
fprintf(stderr, "id = %d, bytes = %lu, size = %d\n", g_proc_id, bytes, 4*sizeof(su3)/8);
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
}
#ifdef MPI
else {
if(g_cart_id == id){
memcpy(&tmp3[0], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][1], sizeof(su3));
memcpy(&tmp3[1], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][2], sizeof(su3));
memcpy(&tmp3[2], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][3], sizeof(su3));
memcpy(&tmp3[3], &g_gauge_field[ g_ipt[tt][X][Y][Z] ][0], sizeof(su3));
# ifndef WORDS_BIGENDIAN
if(prec == 32) {
byte_swap_assign_double2single(tmp2, tmp3, 4*sizeof(su3)/8);
MPI_Send((void*) tmp2, 4*sizeof(su3)/8, MPI_FLOAT, 0, tag, g_cart_grid);
}
else {
byte_swap_assign(tmp, tmp3, 4*sizeof(su3)/8);
MPI_Send((void*) tmp, 4*sizeof(su3)/8, MPI_DOUBLE, 0, tag, g_cart_grid);
}
# else
if(prec == 32) {
double2single(tmp2, tmp3, 4*sizeof(su3)/8);
MPI_Send((void*) tmp2, 4*sizeof(su3)/8, MPI_FLOAT, 0, tag, g_cart_grid);
}
else {
MPI_Send((void*) tmp3, 4*sizeof(su3)/8, MPI_DOUBLE, 0, tag, g_cart_grid);
}
# endif
}
}
#endif
tag++;
}
#ifdef MPI
MPI_Barrier(g_cart_grid);
#endif
}
}
}
return(0);
}
int read_binary_gauge_data(LimeReader * limereader,
const double prec, DML_Checksum * ans) {
int t, x, y , z, status=0;
/* n_uint64_t bytes; */
n_uint64_t bytes;
su3 tmp[4];
float tmp2[72];
DML_SiteRank rank;
DML_checksum_init(ans);
if(prec == 32) bytes = (n_uint64_t) 4*sizeof(su3)/2;
else bytes = (n_uint64_t) 4*sizeof(su3);
for(t = 0; t < T; t++){
for(z = 0; z < LZ; z++){
for(y = 0; y < LY; y++){
#if (defined MPI)
limeReaderSeek(limereader,(n_uint64_t)
(((n_uint64_t) g_proc_coords[1]*LX) +
((n_uint64_t) (((g_proc_coords[0]*T+t)*g_nproc_z*LZ+g_proc_coords[3]*LZ+z)*g_nproc_y*LY
+ g_proc_coords[2]*LY+y)*LX*g_nproc_x))*bytes,
SEEK_SET);
#endif
for(x = 0; x < LX; x++) {
rank = (DML_SiteRank) (g_proc_coords[1]*LX +
(((g_proc_coords[0]*T+t)*g_nproc_z*LZ+g_proc_coords[3]*LZ+z)*g_nproc_y*LY
+ g_proc_coords[2]*LY+y)*((DML_SiteRank)LX*g_nproc_x) + x);
if(prec == 32) {
status = limeReaderReadData(tmp2, &bytes, limereader);
DML_checksum_accum(ans, rank, (char *) tmp2, bytes);
}
else {
status = limeReaderReadData(tmp, &bytes, limereader);
DML_checksum_accum(ans, rank, (char *) tmp, bytes);
}
if(status < 0 && status != LIME_EOR) {
fprintf(stderr, "LIME read error occured with status = %d while reading!\n Aborting...\n", status);
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
#ifndef WORDS_BIGENDIAN
if(prec == 32) {
byte_swap_assign_single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][0], &tmp2[3*18], sizeof(su3)/8);
byte_swap_assign_single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][1], &tmp2[0*18], sizeof(su3)/8);
byte_swap_assign_single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][2], &tmp2[1*18], sizeof(su3)/8);
byte_swap_assign_single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][3], &tmp2[2*18], sizeof(su3)/8);
}
else {
byte_swap_assign(&g_gauge_field[ g_ipt[t][x][y][z] ][0], &tmp[3], sizeof(su3)/8);
byte_swap_assign(&g_gauge_field[ g_ipt[t][x][y][z] ][1], &tmp[0], sizeof(su3)/8);
byte_swap_assign(&g_gauge_field[ g_ipt[t][x][y][z] ][2], &tmp[1], sizeof(su3)/8);
byte_swap_assign(&g_gauge_field[ g_ipt[t][x][y][z] ][3], &tmp[2], sizeof(su3)/8);
}
#else
if(prec == 32) {
single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][0], &tmp2[3*18], sizeof(su3)/8);
single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][1], &tmp2[0], sizeof(su3)/8);
single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][2], &tmp2[18], sizeof(su3)/8);
single2double(&g_gauge_field[ g_ipt[t][x][y][z] ][3], &tmp2[2*18], sizeof(su3)/8);
}
else {
memcpy(&g_gauge_field[ g_ipt[t][x][y][z] ][0], &tmp[3], sizeof(su3));
memcpy(&g_gauge_field[ g_ipt[t][x][y][z] ][1], &tmp[0], sizeof(su3));
memcpy(&g_gauge_field[ g_ipt[t][x][y][z] ][2], &tmp[1], sizeof(su3));
memcpy(&g_gauge_field[ g_ipt[t][x][y][z] ][3], &tmp[2], sizeof(su3));
}
#endif
}
}
}
}
#ifdef MPI
DML_checksum_combine(ans);
#endif
return(0);
}
int write_lime_gauge_field(char * filename, const double plaq, const int counter, const int prec) {
FILE * ofs = NULL;
LimeWriter * limewriter = NULL;
LimeRecordHeader * limeheader = NULL;
/* Message end and Message begin flag */
int ME_flag=0, MB_flag=0, status=0;
n_uint64_t bytes;
DML_Checksum checksum;
write_xlf_info(plaq, counter, filename, 0);
if(g_cart_id == 0) {
ofs = fopen(filename, "a");
if(ofs == (FILE*)NULL) {
fprintf(stderr, "Could not open file %s for writing!\n Aboring...\n", filename);
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
limewriter = limeCreateWriter( ofs );
if(limewriter == (LimeWriter*)NULL) {
fprintf(stderr, "LIME error in file %s for writing!\n Aboring...\n", filename);
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
write_ildg_format_xml("temp.xml", limewriter, prec);
bytes = ((n_uint64_t)LX*g_nproc_x)*((n_uint64_t)LY*g_nproc_y)*((n_uint64_t)LZ*g_nproc_z)*((n_uint64_t)T*g_nproc_t)*((n_uint64_t)4*sizeof(su3));
if(prec == 32) bytes = bytes/((n_uint64_t)2);
MB_flag=0; ME_flag=0;
limeheader = limeCreateHeader(MB_flag, ME_flag, "ildg-binary-data", bytes);
status = limeWriteRecordHeader( limeheader, limewriter);
if(status < 0 ) {
fprintf(stderr, "LIME write header error %d\n", status);
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
limeDestroyHeader( limeheader );
}
write_binary_gauge_data(limewriter, prec, &checksum);
if(g_proc_id == 0 && g_debug_level > 1) {
printf("# checksum for Gauge field written to file %s is %#x %#x\n",
filename, checksum.suma, checksum.sumb);
}
if(g_cart_id == 0) {
limeDestroyWriter( limewriter );
fflush(ofs);
fclose(ofs);
}
write_checksum(filename, &checksum);
return(0);
}
int read_lime_gauge_field(char * filename) {
FILE * ifs;
int status, prec = 64;
n_uint64_t bytes;
char * header_type;
LimeReader * limereader;
DML_Checksum checksum;
ifs = fopen(filename, "r");
if(ifs == (FILE *)NULL) {
fprintf(stderr, "Could not open file %s\n Aborting...\n", filename);
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
limereader = limeCreateReader( ifs );
if( limereader == (LimeReader *)NULL ) {
fprintf(stderr, "Unable to open LimeReader\n");
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(500);
}
while( (status = limeReaderNextRecord(limereader)) != LIME_EOF ) {
if(status != LIME_SUCCESS ) {
fprintf(stderr, "limeReaderNextRecord returned error with status = %d!\n", status);
status = LIME_EOF;
break;
}
header_type = limeReaderType(limereader);
if(strcmp("ildg-binary-data",header_type) == 0) break;
}
if(status == LIME_EOF) {
if(g_proc_id == 0) {
fprintf(stderr, "no ildg-binary-data record found in file %s\n",filename);
fprintf(stderr, "trying old deprecated file format!\n");
}
limeDestroyReader(limereader);
fclose(ifs);
read_gauge_field_time_p(filename);
return(0);
}
bytes = limeReaderBytes(limereader);
if(bytes == ((n_uint64_t)LX*g_nproc_x)*((n_uint64_t)LY*g_nproc_y)*((n_uint64_t)LZ*g_nproc_z)*((n_uint64_t)T*g_nproc_t)*((n_uint64_t)4*sizeof(su3))) prec = 64;
else if(bytes == ((n_uint64_t)LX*g_nproc_x)*((n_uint64_t)LY*g_nproc_y)*((n_uint64_t)LZ*g_nproc_z)*((n_uint64_t)T*g_nproc_t)*((n_uint64_t)4*sizeof(su3)/2)) prec = 32;
else {
fprintf(stderr, "Probably wrong lattice size or precision (bytes=%lu) in file %s\n", bytes, filename);
fprintf(stderr, "Aborting...!\n");
fflush( stdout );
#ifdef MPI
MPI_Abort(MPI_COMM_WORLD, 1);
MPI_Finalize();
#endif
exit(501);
}
if(g_proc_id == 0 && g_debug_level > 2) {
printf("# %d Bit precision read\n", prec);
}
read_binary_gauge_data(limereader, prec, &checksum);
limeDestroyReader(limereader);
fclose(ifs);
if(g_proc_id == 0 && g_debug_level > 1) {
printf("# checksum for gaugefield %s is %#x %#x\n",
filename, checksum.suma, checksum.sumb);
}
return(0);
}